Actual source code: dmpleximpl.h

  1: #pragma once

  3: #include <petscmat.h>
  4: #include <petscdmplex.h>
  5: #include <petscdmplextransform.h>
  6: #include <petscbt.h>
  7: #include <petscsf.h>
  8: #include <petsc/private/dmimpl.h>

 10: #if defined(PETSC_HAVE_EXODUSII)
 11:   #include <exodusII.h>
 12: #endif

 14: PETSC_EXTERN PetscBool  Plexcite;
 15: PETSC_EXTERN const char PlexCitation[];

 17: PETSC_EXTERN PetscLogEvent DMPLEX_Interpolate;
 18: PETSC_EXTERN PetscLogEvent DMPLEX_Partition;
 19: PETSC_EXTERN PetscLogEvent DMPLEX_PartSelf;
 20: PETSC_EXTERN PetscLogEvent DMPLEX_PartLabelInvert;
 21: PETSC_EXTERN PetscLogEvent DMPLEX_PartLabelCreateSF;
 22: PETSC_EXTERN PetscLogEvent DMPLEX_PartStratSF;
 23: PETSC_EXTERN PetscLogEvent DMPLEX_CreatePointSF;
 24: PETSC_EXTERN PetscLogEvent DMPLEX_Distribute;
 25: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeCones;
 26: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeLabels;
 27: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeSF;
 28: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeOverlap;
 29: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeField;
 30: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeData;
 31: PETSC_EXTERN PetscLogEvent DMPLEX_Migrate;
 32: PETSC_EXTERN PetscLogEvent DMPLEX_InterpolateSF;
 33: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalToNaturalBegin;
 34: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalToNaturalEnd;
 35: PETSC_EXTERN PetscLogEvent DMPLEX_NaturalToGlobalBegin;
 36: PETSC_EXTERN PetscLogEvent DMPLEX_NaturalToGlobalEnd;
 37: PETSC_EXTERN PetscLogEvent DMPLEX_Stratify;
 38: PETSC_EXTERN PetscLogEvent DMPLEX_Symmetrize;
 39: PETSC_EXTERN PetscLogEvent DMPLEX_Preallocate;
 40: PETSC_EXTERN PetscLogEvent DMPLEX_ResidualFEM;
 41: PETSC_EXTERN PetscLogEvent DMPLEX_JacobianFEM;
 42: PETSC_EXTERN PetscLogEvent DMPLEX_InterpolatorFEM;
 43: PETSC_EXTERN PetscLogEvent DMPLEX_InjectorFEM;
 44: PETSC_EXTERN PetscLogEvent DMPLEX_IntegralFEM;
 45: PETSC_EXTERN PetscLogEvent DMPLEX_CreateGmsh;
 46: PETSC_EXTERN PetscLogEvent DMPLEX_RebalanceSharedPoints;
 47: PETSC_EXTERN PetscLogEvent DMPLEX_CreateFromFile;
 48: PETSC_EXTERN PetscLogEvent DMPLEX_CreateFromOptions;
 49: PETSC_EXTERN PetscLogEvent DMPLEX_BuildFromCellList;
 50: PETSC_EXTERN PetscLogEvent DMPLEX_BuildCoordinatesFromCellList;
 51: PETSC_EXTERN PetscLogEvent DMPLEX_LocatePoints;
 52: PETSC_EXTERN PetscLogEvent DMPLEX_TopologyView;
 53: PETSC_EXTERN PetscLogEvent DMPLEX_DistributionView;
 54: PETSC_EXTERN PetscLogEvent DMPLEX_LabelsView;
 55: PETSC_EXTERN PetscLogEvent DMPLEX_CoordinatesView;
 56: PETSC_EXTERN PetscLogEvent DMPLEX_SectionView;
 57: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalVectorView;
 58: PETSC_EXTERN PetscLogEvent DMPLEX_LocalVectorView;
 59: PETSC_EXTERN PetscLogEvent DMPLEX_TopologyLoad;
 60: PETSC_EXTERN PetscLogEvent DMPLEX_DistributionLoad;
 61: PETSC_EXTERN PetscLogEvent DMPLEX_LabelsLoad;
 62: PETSC_EXTERN PetscLogEvent DMPLEX_CoordinatesLoad;
 63: PETSC_EXTERN PetscLogEvent DMPLEX_SectionLoad;
 64: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalVectorLoad;
 65: PETSC_EXTERN PetscLogEvent DMPLEX_LocalVectorLoad;
 66: PETSC_EXTERN PetscLogEvent DMPLEX_MetricEnforceSPD;
 67: PETSC_EXTERN PetscLogEvent DMPLEX_MetricNormalize;
 68: PETSC_EXTERN PetscLogEvent DMPLEX_MetricAverage;
 69: PETSC_EXTERN PetscLogEvent DMPLEX_MetricIntersection;
 70: PETSC_EXTERN PetscLogEvent DMPLEX_Generate;
 71: PETSC_EXTERN PetscLogEvent DMPLEX_Transform;
 72: PETSC_EXTERN PetscLogEvent DMPLEX_GetLocalOffsets;

 74: PETSC_EXTERN PetscLogEvent DMPLEX_RebalBuildGraph;
 75: PETSC_EXTERN PetscLogEvent DMPLEX_RebalRewriteSF;
 76: PETSC_EXTERN PetscLogEvent DMPLEX_RebalGatherGraph;
 77: PETSC_EXTERN PetscLogEvent DMPLEX_RebalPartition;
 78: PETSC_EXTERN PetscLogEvent DMPLEX_RebalScatterPart;
 79: PETSC_EXTERN PetscLogEvent DMPLEX_Uninterpolate;

 81: /* Utility struct to store the contents of a Fluent file in memory */
 82: typedef struct {
 83:   int          index; /* Type of section */
 84:   unsigned int zoneID;
 85:   unsigned int first;
 86:   unsigned int last;
 87:   int          type;
 88:   int          nd; /* Either ND or element-type */
 89:   void        *data;
 90: } FluentSection;

 92: struct _PetscGridHash {
 93:   PetscInt     dim;
 94:   PetscReal    lower[3];    /* The lower-left corner */
 95:   PetscReal    upper[3];    /* The upper-right corner */
 96:   PetscReal    extent[3];   /* The box size */
 97:   PetscReal    h[3];        /* The subbox size */
 98:   PetscInt     n[3];        /* The number of subboxes */
 99:   PetscSection cellSection; /* Offsets for cells in each subbox*/
100:   IS           cells;       /* List of cells in each subbox */
101:   DMLabel      cellsSparse; /* Sparse storage for cell map */
102: };

104: typedef struct {
105:   PetscBool isotropic;               /* Is the metric isotropic? */
106:   PetscBool uniform;                 /* Is the metric uniform? */
107:   PetscBool restrictAnisotropyFirst; /* Should anisotropy or normalization come first? */
108:   PetscBool noInsert;                /* Should node insertion/deletion be turned off? */
109:   PetscBool noSwap;                  /* Should facet swapping be turned off? */
110:   PetscBool noMove;                  /* Should node movement be turned off? */
111:   PetscBool noSurf;                  /* Should surface modification be turned off? */
112:   PetscReal h_min, h_max;            /* Minimum/maximum tolerated metric magnitudes */
113:   PetscReal a_max;                   /* Maximum tolerated anisotropy */
114:   PetscReal targetComplexity;        /* Target metric complexity */
115:   PetscReal p;                       /* Degree for L-p normalization methods */
116:   PetscReal gradationFactor;         /* Maximum tolerated length ratio for adjacent edges */
117:   PetscReal hausdorffNumber;         /* Max. distance between piecewise linear representation of boundary and reconstructed ideal boundary */
118:   PetscInt  numIter;                 /* Number of ParMmg mesh adaptation iterations */
119:   PetscInt  verbosity;               /* Level of verbosity for remesher (-1 = no output, 10 = maximum) */
120: } DMPlexMetricCtx;

122: /* Point Numbering in Plex:

124:    Points are numbered contiguously by stratum. Strate are organized as follows:

126:    First Stratum:  Cells [height 0]
127:    Second Stratum: Vertices [depth 0]
128:    Third Stratum:  Faces [height 1]
129:    Fourth Stratum: Edges [depth 1]

131:    We do this so that the numbering of a cell-vertex mesh does not change after interpolation. Within a given stratum,
132:    we allow additional segregation of by cell type.
133: */
134: typedef struct {
135:   PetscInt refct;

137:   PetscSection coneSection;      /* Layout of cones (inedges for DAG) */
138:   PetscInt    *cones;            /* Cone for each point */
139:   PetscInt    *coneOrientations; /* Orientation of each cone point, means cone traversal should start on point 'o', and if negative start on -(o+1) and go in reverse */
140:   PetscSection supportSection;   /* Layout of cones (inedges for DAG) */
141:   PetscInt    *supports;         /* Cone for each point */

143:   struct {                  // DMPolytopeType is an enum (usually size 4), but this needs frequent access
144:     uint8_t value_as_uint8; // in a struct to guard for stronger typing
145:   } *cellTypes;

147:   /* Transformation */
148:   DMPlexTransform tr;                                               /* Type of transform used to define an ephemeral mesh */
149:   char           *transformType;                                    /* Type of transform for uniform cell refinement */
150:   PetscBool       refinementUniform;                                /* Flag for uniform cell refinement */
151:   PetscReal       refinementLimit;                                  /* Maximum volume for refined cell */
152:   PetscErrorCode (*refinementFunc)(const PetscReal[], PetscReal *); /* Function giving the maximum volume for refined cell */

154:   /* Interpolation */
155:   DMPlexInterpolatedFlag interpolated;
156:   DMPlexInterpolatedFlag interpolatedCollective;

158:   /* Ordering */
159:   DMReorderDefaultFlag reorderDefault; /* Reorder the DM by default */

161:   /* Distribution */
162:   PetscBool distDefault;      /* Distribute the DM by default */
163:   PetscInt  overlap;          /* Overlap of the partitions as passed to DMPlexDistribute() or DMPlexDistributeOverlap() */
164:   PetscInt  numOvLabels;      /* The number of labels used for candidate overlap points */
165:   DMLabel   ovLabels[16];     /* Labels used for candidate overlap points */
166:   PetscInt  ovValues[16];     /* Label values used for candidate overlap points */
167:   PetscInt  numOvExLabels;    /* The number of labels used for exclusion */
168:   DMLabel   ovExLabels[16];   /* Labels used to exclude points from the overlap */
169:   PetscInt  ovExValues[16];   /* Label values used to exclude points from the overlap */
170:   char     *distributionName; /* Name of the specific parallel distribution of the DM */
171:   MPI_Comm  nonempty_comm;    /* Communicator used for visualization when some processes do not have cells */

173:   /* Hierarchy */
174:   PetscBool regularRefinement; /* This flag signals that we are a regular refinement of coarseMesh */

176:   /* Generation */
177:   char            *tetgenOpts;
178:   char            *triangleOpts;
179:   PetscPartitioner partitioner;
180:   PetscBool        partitionBalance; /* Evenly divide partition overlap when distributing */
181:   PetscBool        remeshBd;

183:   /* Submesh */
184:   DMLabel          subpointMap;   /* Label each original mesh point in the submesh with its depth, subpoint are the implicit numbering */
185:   IS               subpointIS;    /* IS holding point number in the enclosing mesh of every point in the submesh chart */
186:   PetscObjectState subpointState; /* The state of subpointMap when the subpointIS was last created */

188:   /* Labels and numbering */
189:   PetscObjectState depthState;    /* State of depth label, so that we can determine if a user changes it */
190:   PetscObjectState celltypeState; /* State of celltype label, so that we can determine if a user changes it */
191:   IS               globalVertexNumbers;
192:   IS               globalCellNumbers;

194:   /* Constraints */
195:   PetscSection anchorSection;          /* maps constrained points to anchor points */
196:   IS           anchorIS;               /* anchors indexed by the above section */
197:   PetscErrorCode (*createanchors)(DM); /* automatically compute anchors (probably from tree constraints) */
198:   PetscErrorCode (*computeanchormatrix)(DM, PetscSection, PetscSection, Mat);

200:   /* Tree: automatically construct constraints for hierarchically non-conforming meshes */
201:   PetscSection parentSection; /* dof == 1 if point has parent */
202:   PetscInt    *parents;       /* point to parent */
203:   PetscInt    *childIDs;      /* point to child ID */
204:   PetscSection childSection;  /* inverse of parent section */
205:   PetscInt    *children;      /* point to children */
206:   DM           referenceTree; /* reference tree to which child ID's refer */
207:   PetscErrorCode (*getchildsymmetry)(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *);

209:   /* MATIS support */
210:   PetscSection subdomainSection;

212:   /* Adjacency */
213:   PetscBool useAnchors;                                                          /* Replace constrained points with their anchors in adjacency lists */
214:   PetscErrorCode (*useradjacency)(DM, PetscInt, PetscInt *, PetscInt[], void *); /* User callback for adjacency */
215:   void *useradjacencyctx;                                                        /* User context for callback */

217:   // Periodicity
218:   struct {
219:     // Specified by the user
220:     PetscInt num_face_sfs;          // number of face_sfs
221:     PetscSF *face_sfs;              // root(donor faces) <-- leaf(local faces)
222:     PetscScalar (*transform)[4][4]; // geometric transform
223:     // Created eagerly (depends on points)
224:     PetscSF composed_sf; // root(non-periodic global points) <-- leaf(local points)
225:     IS     *periodic_points;
226:   } periodic;

228:   /* Projection */
229:   PetscInt maxProjectionHeight; /* maximum height of cells used in DMPlexProject functions */
230:   PetscInt activePoint;         /* current active point in iteration */

232:   /* Output */
233:   PetscInt  vtkCellHeight;          /* The height of cells for output, default is 0 */
234:   PetscReal scale[NUM_PETSC_UNITS]; /* The scale for each SI unit */

236:   /* Geometry */
237:   PetscReal     minradius;                        /* Minimum distance from cell centroid to face */
238:   PetscBool     useHashLocation;                  /* Use grid hashing for point location */
239:   PetscGridHash lbox;                             /* Local box for searching */
240:   void (*coordFunc)(PetscInt, PetscInt, PetscInt, /* Function used to remap newly introduced vertices */
241:                     const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);

243:   /* Neighbors */
244:   PetscMPIInt *neighbors;

246:   /* Metric */
247:   DMPlexMetricCtx *metricCtx;

249:   /* FEM */
250:   PetscBool useCeed;      /* This should convert to a registration system when there are more FEM backends */
251:   PetscBool useMatClPerm; /* Use the closure permutation when assembling matrices */

253:   /* Debugging */
254:   PetscBool printSetValues;
255:   PetscInt  printFEM;
256:   PetscInt  printFVM;
257:   PetscInt  printL2;
258:   PetscInt  printLocate;
259:   PetscReal printTol;
260: } DM_Plex;

262: PETSC_INTERN PetscErrorCode DMPlexCopy_Internal(DM, PetscBool, PetscBool, DM);
263: PETSC_INTERN PetscErrorCode DMPlexReplace_Internal(DM, DM *);

265: PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll_VTU(DM, PetscViewer);
266: PETSC_EXTERN PetscErrorCode VecView_Plex_Local(Vec, PetscViewer);
267: PETSC_EXTERN PetscErrorCode VecView_Plex_Native(Vec, PetscViewer);
268: PETSC_EXTERN PetscErrorCode VecView_Plex(Vec, PetscViewer);
269: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Local(Vec, PetscViewer);
270: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Native(Vec, PetscViewer);
271: PETSC_EXTERN PetscErrorCode VecLoad_Plex(Vec, PetscViewer);
272: PETSC_INTERN PetscErrorCode DMPlexGetFieldType_Internal(DM, PetscSection, PetscInt, PetscInt *, PetscInt *, PetscViewerVTKFieldType *);
273: PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM, PetscViewer);
274: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject, PetscViewer);
275: #if defined(PETSC_HAVE_HDF5)
276: PETSC_EXTERN PetscErrorCode VecView_Plex_Local_HDF5(Vec, PetscViewer);
277: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5(Vec, PetscViewer);
278: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5(Vec, PetscViewer);
279: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5_Native(Vec, PetscViewer);
280: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5_Native(Vec, PetscViewer);
281: PETSC_EXTERN PetscErrorCode DMPlexView_HDF5(DM, PetscViewer);
282: PETSC_EXTERN PetscErrorCode DMPlexLoad_HDF5(DM, PetscViewer);
283: PETSC_INTERN PetscErrorCode DMPlexTopologyView_HDF5_Internal(DM, IS, PetscViewer);
284: PETSC_INTERN PetscErrorCode DMPlexCoordinatesView_HDF5_Internal(DM, PetscViewer);
285: PETSC_INTERN PetscErrorCode DMPlexLabelsView_HDF5_Internal(DM, IS, PetscViewer);
286: PETSC_INTERN PetscErrorCode DMPlexSectionView_HDF5_Internal(DM, PetscViewer, DM);
287: PETSC_INTERN PetscErrorCode DMPlexGlobalVectorView_HDF5_Internal(DM, PetscViewer, DM, Vec);
288: PETSC_INTERN PetscErrorCode DMPlexLocalVectorView_HDF5_Internal(DM, PetscViewer, DM, Vec);
289: PETSC_INTERN PetscErrorCode DMPlexTopologyLoad_HDF5_Internal(DM, PetscViewer, PetscSF *);
290: PETSC_INTERN PetscErrorCode DMPlexCoordinatesLoad_HDF5_Internal(DM, PetscViewer, PetscSF);
291: PETSC_INTERN PetscErrorCode DMPlexLabelsLoad_HDF5_Internal(DM, PetscViewer, PetscSF);
292: PETSC_INTERN PetscErrorCode DMPlexSectionLoad_HDF5_Internal(DM, PetscViewer, DM, PetscSF, PetscSF *, PetscSF *);
293: PETSC_INTERN PetscErrorCode DMPlexVecLoad_HDF5_Internal(DM, PetscViewer, DM, PetscSF, Vec);
294: PETSC_INTERN PetscErrorCode DMPlexView_HDF5_Internal(DM, PetscViewer);
295: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Internal(DM, PetscViewer);
296: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Xdmf_Internal(DM, PetscViewer);
297: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Internal(Vec, PetscViewer);
298: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Native_Internal(Vec, PetscViewer);
299: PETSC_INTERN PetscErrorCode VecView_Plex_Local_HDF5_Internal(Vec, PetscViewer);
300: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Internal(Vec, PetscViewer);
301: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Native_Internal(Vec, PetscViewer);
302: #endif
303: PETSC_EXTERN PetscErrorCode VecView_Plex_Local_CGNS(Vec, PetscViewer);
304: #if defined(PETSC_HAVE_CGNS)
305: PETSC_EXTERN PetscErrorCode VecLoad_Plex_CGNS_Internal(Vec, PetscViewer);
306: #endif

308: PETSC_INTERN PetscErrorCode DMPlexVecGetClosure_Internal(DM, PetscSection, PetscBool, Vec, PetscInt, PetscInt *, PetscScalar *[]);
309: PETSC_INTERN PetscErrorCode DMPlexVecGetClosureAtDepth_Internal(DM, PetscSection, Vec, PetscInt, PetscInt, PetscInt *, PetscScalar *[]);
310: PETSC_INTERN PetscErrorCode DMPlexClosurePoints_Private(DM, PetscInt, const PetscInt[], IS *);
311: PETSC_INTERN PetscErrorCode DMSetFromOptions_NonRefinement_Plex(DM, PetscOptionItems *);
312: PETSC_INTERN PetscErrorCode DMSetFromOptions_Overlap_Plex(DM, PetscOptionItems *, PetscInt *);
313: PETSC_INTERN PetscErrorCode DMCoarsen_Plex(DM, MPI_Comm, DM *);
314: PETSC_INTERN PetscErrorCode DMCoarsenHierarchy_Plex(DM, PetscInt, DM[]);
315: PETSC_INTERN PetscErrorCode DMRefine_Plex(DM, MPI_Comm, DM *);
316: PETSC_INTERN PetscErrorCode DMRefineHierarchy_Plex(DM, PetscInt, DM[]);
317: PETSC_INTERN PetscErrorCode DMAdaptLabel_Plex(DM, Vec, DMLabel, DMLabel, DM *);
318: PETSC_INTERN PetscErrorCode DMExtrude_Plex(DM, PetscInt, DM *);
319: PETSC_INTERN PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
320: PETSC_INTERN PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
321: PETSC_INTERN PetscErrorCode DMProjectFunctionLocal_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
322: PETSC_INTERN PetscErrorCode DMProjectFunctionLabelLocal_Plex(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
323: PETSC_INTERN PetscErrorCode DMProjectFieldLocal_Plex(DM, PetscReal, Vec, void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode, Vec);
324: PETSC_INTERN PetscErrorCode DMProjectFieldLabelLocal_Plex(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Vec, void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode, Vec);
325: PETSC_INTERN PetscErrorCode DMProjectBdFieldLabelLocal_Plex(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Vec, void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode, Vec);
326: PETSC_INTERN PetscErrorCode DMComputeL2Diff_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
327: PETSC_INTERN PetscErrorCode DMComputeL2GradientDiff_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal[], PetscReal *);
328: PETSC_INTERN PetscErrorCode DMComputeL2FieldDiff_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
329: PETSC_INTERN PetscErrorCode DMLocatePoints_Plex(DM, Vec, DMPointLocationType, PetscSF);

331: #if defined(PETSC_HAVE_EXODUSII)
332: PETSC_INTERN PetscErrorCode DMView_PlexExodusII(DM, PetscViewer);
333: PETSC_INTERN PetscErrorCode VecView_PlexExodusII_Internal(Vec, PetscViewer);
334: PETSC_INTERN PetscErrorCode VecLoad_PlexExodusII_Internal(Vec, PetscViewer);
335: #endif
336: PETSC_INTERN PetscErrorCode DMView_PlexCGNS(DM, PetscViewer);
337: PETSC_INTERN PetscErrorCode DMPlexCreateCGNSFromFile_Internal(MPI_Comm, const char *, PetscBool, DM *);
338: PETSC_INTERN PetscErrorCode DMPlexCreateCGNS_Internal_Serial(MPI_Comm, PetscInt, PetscBool, DM *);
339: PETSC_INTERN PetscErrorCode DMPlexCreateCGNS_Internal_Parallel(MPI_Comm, PetscInt, PetscBool, DM *);
340: PETSC_INTERN PetscErrorCode DMPlexVTKGetCellType_Internal(DM, PetscInt, PetscInt, PetscInt *);
341: PETSC_INTERN PetscErrorCode DMPlexGetAdjacency_Internal(DM, PetscInt, PetscBool, PetscBool, PetscBool, PetscInt *, PetscInt *[]);
342: PETSC_INTERN PetscErrorCode DMPlexGetRawFaces_Internal(DM, DMPolytopeType, const PetscInt[], PetscInt *, const DMPolytopeType *[], const PetscInt *[], const PetscInt *[]);
343: PETSC_INTERN PetscErrorCode DMPlexRestoreRawFaces_Internal(DM, DMPolytopeType, const PetscInt[], PetscInt *, const DMPolytopeType *[], const PetscInt *[], const PetscInt *[]);
344: PETSC_INTERN PetscErrorCode DMPlexComputeCellType_Internal(DM, PetscInt, PetscInt, DMPolytopeType *);
345: PETSC_INTERN PetscErrorCode DMPlexVecSetFieldClosure_Internal(DM, PetscSection, Vec, PetscBool[], PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscScalar[], InsertMode);
346: PETSC_INTERN PetscErrorCode DMPlexProjectConstraints_Internal(DM, Vec, Vec);
347: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_SetTree(DM, PetscSection, PetscInt[], PetscInt[]);
348: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_Union(DM, DM, const char *, DM *);
349: PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorTree(DM, DM, PetscSF, PetscInt *, Mat);
350: PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorTree(DM, DM, PetscSF, PetscInt *, Mat);
351: PETSC_INTERN PetscErrorCode DMPlexAnchorsModifyMat(DM, PetscSection, PetscInt, PetscInt, const PetscInt[], const PetscInt ***, const PetscScalar[], PetscInt *, PetscInt *, PetscInt *[], PetscScalar *[], PetscInt[], PetscBool);
352: PETSC_INTERN PetscErrorCode DMPlexAnchorsModifyMat_Internal(DM, PetscSection, PetscInt, PetscInt, const PetscInt[], const PetscInt ***, PetscInt, PetscInt, const PetscScalar[], PetscInt *, PetscInt *, PetscInt *[], PetscScalar *[], PetscInt[], PetscBool, PetscBool);
353: PETSC_INTERN PetscErrorCode DMPlexAnchorsGetSubMatModification(DM, PetscSection, PetscInt, PetscInt, const PetscInt[], const PetscInt ***, PetscInt *, PetscInt *, PetscInt *[], PetscInt[], PetscScalar *[]);
354: PETSC_INTERN PetscErrorCode DMPlexLocatePoint_Internal(DM, PetscInt, const PetscScalar[], PetscInt, PetscInt *);
355: /* this is PETSC_EXTERN just because of src/dm/impls/plex/tests/ex18.c */
356: PETSC_EXTERN PetscErrorCode DMPlexOrientInterface_Internal(DM);
357: PETSC_INTERN PetscErrorCode DMPlexOrientCells_Internal(DM, IS, IS);

359: /* Applications may use this function */
360: PETSC_EXTERN PetscErrorCode DMPlexCreateNumbering_Plex(DM, PetscInt, PetscInt, PetscInt, PetscInt *, PetscSF, IS *);

362: PETSC_INTERN PetscErrorCode DMPlexInterpolateInPlace_Internal(DM);
363: PETSC_INTERN PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Internal(DM, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool);
364: PETSC_INTERN PetscErrorCode DMPlexMigrateIsoperiodicFaceSF_Internal(DM, DM, PetscSF);
365: PETSC_INTERN PetscErrorCode DMPlexCreateVertexNumbering_Internal(DM, PetscBool, IS *);
366: PETSC_INTERN PetscErrorCode DMPlexRefine_Internal(DM, Vec, DMLabel, DMLabel, DM *);
367: PETSC_INTERN PetscErrorCode DMPlexCoarsen_Internal(DM, Vec, DMLabel, DMLabel, DM *);
368: PETSC_INTERN PetscErrorCode DMCreateMatrix_Plex(DM, Mat *);

370: PETSC_INTERN PetscErrorCode DMPlexGetOverlap_Plex(DM, PetscInt *);
371: PETSC_INTERN PetscErrorCode DMPlexSetOverlap_Plex(DM, DM, PetscInt);
372: PETSC_INTERN PetscErrorCode DMPlexDistributeGetDefault_Plex(DM, PetscBool *);
373: PETSC_INTERN PetscErrorCode DMPlexDistributeSetDefault_Plex(DM, PetscBool);
374: PETSC_INTERN PetscErrorCode DMPlexReorderGetDefault_Plex(DM, DMReorderDefaultFlag *);
375: PETSC_INTERN PetscErrorCode DMPlexReorderSetDefault_Plex(DM, DMReorderDefaultFlag);
376: PETSC_INTERN PetscErrorCode DMPlexGetUseCeed_Plex(DM, PetscBool *);
377: PETSC_INTERN PetscErrorCode DMPlexSetUseCeed_Plex(DM, PetscBool);
378: PETSC_INTERN PetscErrorCode DMReorderSectionGetDefault_Plex(DM, DMReorderDefaultFlag *);
379: PETSC_INTERN PetscErrorCode DMReorderSectionSetDefault_Plex(DM, DMReorderDefaultFlag);
380: PETSC_INTERN PetscErrorCode DMReorderSectionGetType_Plex(DM, MatOrderingType *);
381: PETSC_INTERN PetscErrorCode DMReorderSectionSetType_Plex(DM, MatOrderingType);

383: #if 1
384: static inline PetscInt DihedralInvert(PetscInt N, PetscInt a)
385: {
386:   return (a <= 0) ? a : (N - a);
387: }

389: static inline PetscInt DihedralCompose(PetscInt N, PetscInt a, PetscInt b)
390: {
391:   if (!N) return 0;
392:   return (a >= 0) ? ((b >= 0) ? ((a + b) % N) : -(((a - b - 1) % N) + 1)) : ((b >= 0) ? -(((N - b - a - 1) % N) + 1) : ((N + b - a) % N));
393: }

395: static inline PetscInt DihedralSwap(PetscInt N, PetscInt a, PetscInt b)
396: {
397:   return DihedralCompose(N, DihedralInvert(N, a), b);
398: }
399: #else
400: /* TODO
401:    This is a reimplementation of the tensor dihedral symmetries using the new orientations.
402:    These should be turned on when we convert to new-style orientations in p4est.
403: */
404: /* invert dihedral symmetry: return a^-1,
405:  * using the representation described in
406:  * DMPlexGetConeOrientation() */
407: static inline PetscInt DihedralInvert(PetscInt N, PetscInt a)
408: {
409:   switch (N) {
410:   case 0:
411:     return 0;
412:   case 2:
413:     return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_SEGMENT, 0, a);
414:   case 4:
415:     return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_QUADRILATERAL, 0, a);
416:   case 8:
417:     return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_HEXAHEDRON, 0, a);
418:   default:
419:     SETERRABORT(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid celltype for DihedralInvert()");
420:   }
421:   return 0;
422: }

424: /* compose dihedral symmetry: return b * a,
425:  * using the representation described in
426:  * DMPlexGetConeOrientation() */
427: static inline PetscInt DihedralCompose(PetscInt N, PetscInt a, PetscInt b)
428: {
429:   switch (N) {
430:   case 0:
431:     return 0;
432:   case 2:
433:     return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_SEGMENT, b, a);
434:   case 4:
435:     return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_QUADRILATERAL, b, a);
436:   case 8:
437:     return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_HEXAHEDRON, b, a);
438:   default:
439:     SETERRABORT(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid celltype for DihedralCompose()");
440:   }
441:   return 0;
442: }

444: /* swap dihedral symmetries: return b * a^-1,
445:  * using the representation described in
446:  * DMPlexGetConeOrientation() */
447: static inline PetscInt DihedralSwap(PetscInt N, PetscInt a, PetscInt b)
448: {
449:   switch (N) {
450:   case 0:
451:     return 0;
452:   case 2:
453:     return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_SEGMENT, b, a);
454:   case 4:
455:     return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_QUADRILATERAL, b, a);
456:   case 8:
457:     return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_HEXAHEDRON, b, a);
458:   default:
459:     SETERRABORT(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid celltype for DihedralCompose()");
460:   }
461:   return 0;
462: }
463: #endif
464: PETSC_EXTERN PetscInt       DMPolytopeConvertNewOrientation_Internal(DMPolytopeType, PetscInt);
465: PETSC_EXTERN PetscInt       DMPolytopeConvertOldOrientation_Internal(DMPolytopeType, PetscInt);
466: PETSC_EXTERN PetscErrorCode DMPlexConvertOldOrientations_Internal(DM);

468: PETSC_EXTERN PetscErrorCode DMPlexComputeIntegral_Internal(DM, Vec, PetscInt, PetscInt, PetscScalar *, void *);
469: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Internal(DM, PetscFormKey, IS, PetscReal, Vec, Vec, PetscReal, Vec, void *);
470: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Hybrid_Internal(DM, PetscFormKey[], IS, PetscReal, Vec, Vec, PetscReal, Vec, void *);
471: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Internal(DM, PetscFormKey, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
472: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Hybrid_Internal(DM, PetscFormKey[], IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
473: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Action_Internal(DM, PetscFormKey, IS, PetscReal, PetscReal, Vec, Vec, Vec, Vec, void *);
474: PETSC_EXTERN PetscErrorCode DMPlexReconstructGradients_Internal(DM, PetscFV, PetscInt, PetscInt, Vec, Vec, Vec, Vec);

476: /* Matvec with A in row-major storage, x and y can be aliased */
477: static inline void DMPlex_Mult2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
478: {
479:   const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
480:   y[0 * ldx]             = A[0] * z[0] + A[1] * z[1];
481:   y[1 * ldx]             = A[2] * z[0] + A[3] * z[1];
482:   (void)PetscLogFlops(6.0);
483: }
484: static inline void DMPlex_Mult3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
485: {
486:   const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
487:   y[0 * ldx]             = A[0] * z[0] + A[1] * z[1] + A[2] * z[2];
488:   y[1 * ldx]             = A[3] * z[0] + A[4] * z[1] + A[5] * z[2];
489:   y[2 * ldx]             = A[6] * z[0] + A[7] * z[1] + A[8] * z[2];
490:   (void)PetscLogFlops(15.0);
491: }
492: static inline void DMPlex_MultTranspose2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
493: {
494:   const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
495:   y[0 * ldx]             = A[0] * z[0] + A[2] * z[1];
496:   y[1 * ldx]             = A[1] * z[0] + A[3] * z[1];
497:   (void)PetscLogFlops(6.0);
498: }
499: static inline void DMPlex_MultTranspose3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
500: {
501:   const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
502:   y[0 * ldx]             = PetscConj(A[0]) * z[0] + PetscConj(A[3]) * z[1] + PetscConj(A[6]) * z[2];
503:   y[1 * ldx]             = PetscConj(A[1]) * z[0] + PetscConj(A[4]) * z[1] + PetscConj(A[7]) * z[2];
504:   y[2 * ldx]             = PetscConj(A[2]) * z[0] + PetscConj(A[5]) * z[1] + PetscConj(A[8]) * z[2];
505:   (void)PetscLogFlops(15.0);
506: }
507: static inline void DMPlex_Mult2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
508: {
509:   const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
510:   y[0 * ldx]             = A[0] * z[0] + A[1] * z[1];
511:   y[1 * ldx]             = A[2] * z[0] + A[3] * z[1];
512:   (void)PetscLogFlops(6.0);
513: }
514: static inline void DMPlex_Mult3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
515: {
516:   const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
517:   y[0 * ldx]             = A[0] * z[0] + A[1] * z[1] + A[2] * z[2];
518:   y[1 * ldx]             = A[3] * z[0] + A[4] * z[1] + A[5] * z[2];
519:   y[2 * ldx]             = A[6] * z[0] + A[7] * z[1] + A[8] * z[2];
520:   (void)PetscLogFlops(15.0);
521: }
522: static inline void DMPlex_MultAdd2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
523: {
524:   const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
525:   y[0 * ldx] += A[0] * z[0] + A[1] * z[1];
526:   y[1 * ldx] += A[2] * z[0] + A[3] * z[1];
527:   (void)PetscLogFlops(6.0);
528: }
529: static inline void DMPlex_MultAdd3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
530: {
531:   const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
532:   y[0 * ldx] += A[0] * z[0] + A[1] * z[1] + A[2] * z[2];
533:   y[1 * ldx] += A[3] * z[0] + A[4] * z[1] + A[5] * z[2];
534:   y[2 * ldx] += A[6] * z[0] + A[7] * z[1] + A[8] * z[2];
535:   (void)PetscLogFlops(15.0);
536: }
537: /*
538:   A: packed, row-major m x n array
539:   x: length m array
540:   y: length n arra
541:   ldx: the stride in x and y

543:   A[i,j] = A[i * n + j]
544:   A^T[j,i] = A[i,j]
545: */
546: static inline void DMPlex_MultTransposeReal_Internal(const PetscReal A[], PetscInt m, PetscInt n, PetscInt ldx, const PetscScalar x[], PetscScalar y[])
547: {
548:   PetscScalar z[3];
549:   PetscInt    i, j;
550:   for (i = 0; i < m; ++i) z[i] = x[i * ldx];
551:   for (j = 0; j < n; ++j) {
552:     const PetscInt l = j * ldx;
553:     y[l]             = 0;
554:     for (i = 0; i < m; ++i) y[l] += A[i * n + j] * z[i];
555:   }
556:   (void)PetscLogFlops(2 * m * n);
557: }
558: static inline void DMPlex_MultTranspose2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
559: {
560:   const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
561:   y[0 * ldx]             = A[0] * z[0] + A[2] * z[1];
562:   y[1 * ldx]             = A[1] * z[0] + A[3] * z[1];
563:   (void)PetscLogFlops(6.0);
564: }
565: static inline void DMPlex_MultTranspose3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
566: {
567:   const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
568:   y[0 * ldx]             = A[0] * z[0] + A[3] * z[1] + A[6] * z[2];
569:   y[1 * ldx]             = A[1] * z[0] + A[4] * z[1] + A[7] * z[2];
570:   y[2 * ldx]             = A[2] * z[0] + A[5] * z[1] + A[8] * z[2];
571:   (void)PetscLogFlops(15.0);
572: }

574: static inline void DMPlex_MatMult2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
575: {
576: #define PLEX_DIM__ 2
577:   PetscScalar z[PLEX_DIM__];
578:   for (PetscInt j = 0; j < n; ++j) {
579:     for (int d = 0; d < PLEX_DIM__; ++d) z[d] = B[d * ldb + j];
580:     DMPlex_Mult2D_Internal(A, 1, z, z);
581:     for (int d = 0; d < PLEX_DIM__; ++d) C[d * ldb + j] = z[d];
582:   }
583:   (void)PetscLogFlops(8.0 * n);
584: #undef PLEX_DIM__
585: }
586: static inline void DMPlex_MatMult3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
587: {
588: #define PLEX_DIM__ 3
589:   PetscScalar z[PLEX_DIM__];
590:   for (PetscInt j = 0; j < n; ++j) {
591:     for (int d = 0; d < PLEX_DIM__; ++d) z[d] = B[d * ldb + j];
592:     DMPlex_Mult3D_Internal(A, 1, z, z);
593:     for (int d = 0; d < PLEX_DIM__; ++d) C[d * ldb + j] = z[d];
594:   }
595:   (void)PetscLogFlops(8.0 * n);
596: #undef PLEX_DIM__
597: }
598: static inline void DMPlex_MatMultTranspose2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
599: {
600: #define PLEX_DIM__ 2
601:   PetscScalar z[PLEX_DIM__];
602:   for (PetscInt j = 0; j < n; ++j) {
603:     for (int d = 0; d < PLEX_DIM__; ++d) z[d] = B[d * ldb + j];
604:     DMPlex_MultTranspose2D_Internal(A, 1, z, z);
605:     for (int d = 0; d < PLEX_DIM__; ++d) C[d * ldb + j] = z[d];
606:   }
607:   (void)PetscLogFlops(8.0 * n);
608: #undef PLEX_DIM__
609: }
610: static inline void DMPlex_MatMultTranspose3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
611: {
612: #define PLEX_DIM__ 3
613:   PetscScalar z[PLEX_DIM__];
614:   for (PetscInt j = 0; j < n; ++j) {
615:     for (int d = 0; d < PLEX_DIM__; ++d) z[d] = B[d * ldb + j];
616:     DMPlex_MultTranspose3D_Internal(A, 1, z, z);
617:     for (int d = 0; d < PLEX_DIM__; ++d) C[d * ldb + j] = z[d];
618:   }
619:   (void)PetscLogFlops(8.0 * n);
620: #undef PLEX_DIM__
621: }

623: static inline void DMPlex_MatMultLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
624: {
625:   for (PetscInt j = 0; j < m; ++j) DMPlex_MultTranspose2D_Internal(A, 1, &B[j * ldb], &C[j * ldb]);
626:   (void)PetscLogFlops(8.0 * m);
627: }
628: static inline void DMPlex_MatMultLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
629: {
630:   for (PetscInt j = 0; j < m; ++j) DMPlex_MultTranspose3D_Internal(A, 1, &B[j * ldb], &C[j * ldb]);
631:   (void)PetscLogFlops(8.0 * m);
632: }
633: static inline void DMPlex_MatMultTransposeLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
634: {
635:   for (PetscInt j = 0; j < m; ++j) DMPlex_Mult2D_Internal(A, 1, &B[j * ldb], &C[j * ldb]);
636:   (void)PetscLogFlops(8.0 * m);
637: }
638: static inline void DMPlex_MatMultTransposeLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
639: {
640:   for (PetscInt j = 0; j < m; ++j) DMPlex_Mult3D_Internal(A, 1, &B[j * ldb], &C[j * ldb]);
641:   (void)PetscLogFlops(8.0 * m);
642: }

644: static inline void DMPlex_PTAP2DReal_Internal(const PetscReal P[], const PetscScalar A[], PetscScalar C[])
645: {
646:   PetscScalar out[4];
647:   PetscInt    i, j, k, l;
648:   for (i = 0; i < 2; ++i) {
649:     for (j = 0; j < 2; ++j) {
650:       out[i * 2 + j] = 0.;
651:       for (k = 0; k < 2; ++k) {
652:         for (l = 0; l < 2; ++l) out[i * 2 + j] += P[k * 2 + i] * A[k * 2 + l] * P[l * 2 + j];
653:       }
654:     }
655:   }
656:   for (i = 0; i < 2 * 2; ++i) C[i] = out[i];
657:   (void)PetscLogFlops(48.0);
658: }
659: static inline void DMPlex_PTAP3DReal_Internal(const PetscReal P[], const PetscScalar A[], PetscScalar C[])
660: {
661:   PetscScalar out[9];
662:   PetscInt    i, j, k, l;
663:   for (i = 0; i < 3; ++i) {
664:     for (j = 0; j < 3; ++j) {
665:       out[i * 3 + j] = 0.;
666:       for (k = 0; k < 3; ++k) {
667:         for (l = 0; l < 3; ++l) out[i * 3 + j] += P[k * 3 + i] * A[k * 3 + l] * P[l * 3 + j];
668:       }
669:     }
670:   }
671:   for (i = 0; i < 2 * 2; ++i) C[i] = out[i];
672:   (void)PetscLogFlops(243.0);
673: }
674: /* TODO Fix for aliasing of A and C */
675: static inline void DMPlex_PTAPReal_Internal(const PetscReal P[], PetscInt m, PetscInt n, const PetscScalar A[], PetscScalar C[])
676: {
677:   PetscInt i, j, k, l;
678:   for (i = 0; i < n; ++i) {
679:     for (j = 0; j < n; ++j) {
680:       C[i * n + j] = 0.;
681:       for (k = 0; k < m; ++k) {
682:         for (l = 0; l < m; ++l) C[i * n + j] += P[k * n + i] * A[k * m + l] * P[l * n + j];
683:       }
684:     }
685:   }
686:   (void)PetscLogFlops(243.0);
687: }

689: static inline void DMPlex_Transpose2D_Internal(PetscScalar A[])
690: {
691:   PetscScalar tmp;
692:   tmp  = A[1];
693:   A[1] = A[2];
694:   A[2] = tmp;
695: }
696: static inline void DMPlex_Transpose3D_Internal(PetscScalar A[])
697: {
698:   PetscScalar tmp;
699:   tmp  = A[1];
700:   A[1] = A[3];
701:   A[3] = tmp;
702:   tmp  = A[2];
703:   A[2] = A[6];
704:   A[6] = tmp;
705:   tmp  = A[5];
706:   A[5] = A[7];
707:   A[7] = tmp;
708: }

710: static inline void DMPlex_Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
711: {
712:   // Allow zero volume cells
713:   const PetscReal invDet = detJ == 0 ? 1.0 : (PetscReal)1.0 / detJ;

715:   invJ[0] = invDet * J[3];
716:   invJ[1] = -invDet * J[1];
717:   invJ[2] = -invDet * J[2];
718:   invJ[3] = invDet * J[0];
719:   (void)PetscLogFlops(5.0);
720: }

722: static inline void DMPlex_Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
723: {
724:   // Allow zero volume cells
725:   const PetscReal invDet = detJ == 0 ? 1.0 : (PetscReal)1.0 / detJ;

727:   invJ[0 * 3 + 0] = invDet * (J[1 * 3 + 1] * J[2 * 3 + 2] - J[1 * 3 + 2] * J[2 * 3 + 1]);
728:   invJ[0 * 3 + 1] = invDet * (J[0 * 3 + 2] * J[2 * 3 + 1] - J[0 * 3 + 1] * J[2 * 3 + 2]);
729:   invJ[0 * 3 + 2] = invDet * (J[0 * 3 + 1] * J[1 * 3 + 2] - J[0 * 3 + 2] * J[1 * 3 + 1]);
730:   invJ[1 * 3 + 0] = invDet * (J[1 * 3 + 2] * J[2 * 3 + 0] - J[1 * 3 + 0] * J[2 * 3 + 2]);
731:   invJ[1 * 3 + 1] = invDet * (J[0 * 3 + 0] * J[2 * 3 + 2] - J[0 * 3 + 2] * J[2 * 3 + 0]);
732:   invJ[1 * 3 + 2] = invDet * (J[0 * 3 + 2] * J[1 * 3 + 0] - J[0 * 3 + 0] * J[1 * 3 + 2]);
733:   invJ[2 * 3 + 0] = invDet * (J[1 * 3 + 0] * J[2 * 3 + 1] - J[1 * 3 + 1] * J[2 * 3 + 0]);
734:   invJ[2 * 3 + 1] = invDet * (J[0 * 3 + 1] * J[2 * 3 + 0] - J[0 * 3 + 0] * J[2 * 3 + 1]);
735:   invJ[2 * 3 + 2] = invDet * (J[0 * 3 + 0] * J[1 * 3 + 1] - J[0 * 3 + 1] * J[1 * 3 + 0]);
736:   (void)PetscLogFlops(37.0);
737: }

739: static inline void DMPlex_Det2D_Internal(PetscReal *detJ, const PetscReal J[])
740: {
741:   *detJ = J[0] * J[3] - J[1] * J[2];
742:   (void)PetscLogFlops(3.0);
743: }

745: static inline void DMPlex_Det3D_Internal(PetscReal *detJ, const PetscReal J[])
746: {
747:   *detJ = (J[0 * 3 + 0] * (J[1 * 3 + 1] * J[2 * 3 + 2] - J[1 * 3 + 2] * J[2 * 3 + 1]) + J[0 * 3 + 1] * (J[1 * 3 + 2] * J[2 * 3 + 0] - J[1 * 3 + 0] * J[2 * 3 + 2]) + J[0 * 3 + 2] * (J[1 * 3 + 0] * J[2 * 3 + 1] - J[1 * 3 + 1] * J[2 * 3 + 0]));
748:   (void)PetscLogFlops(12.0);
749: }

751: static inline void DMPlex_Det2D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
752: {
753:   *detJ = PetscRealPart(J[0]) * PetscRealPart(J[3]) - PetscRealPart(J[1]) * PetscRealPart(J[2]);
754:   (void)PetscLogFlops(3.0);
755: }

757: static inline void DMPlex_Det3D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
758: {
759:   *detJ = (PetscRealPart(J[0 * 3 + 0]) * (PetscRealPart(J[1 * 3 + 1]) * PetscRealPart(J[2 * 3 + 2]) - PetscRealPart(J[1 * 3 + 2]) * PetscRealPart(J[2 * 3 + 1])) + PetscRealPart(J[0 * 3 + 1]) * (PetscRealPart(J[1 * 3 + 2]) * PetscRealPart(J[2 * 3 + 0]) - PetscRealPart(J[1 * 3 + 0]) * PetscRealPart(J[2 * 3 + 2])) + PetscRealPart(J[0 * 3 + 2]) * (PetscRealPart(J[1 * 3 + 0]) * PetscRealPart(J[2 * 3 + 1]) - PetscRealPart(J[1 * 3 + 1]) * PetscRealPart(J[2 * 3 + 0])));
760:   (void)PetscLogFlops(12.0);
761: }

763: static inline void DMPlex_WaxpyD_Internal(PetscInt dim, PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w)
764: {
765:   PetscInt d;
766:   for (d = 0; d < dim; ++d) w[d] = a * x[d] + y[d];
767: }

769: static inline PetscReal DMPlex_DotD_Internal(PetscInt dim, const PetscScalar *x, const PetscReal *y)
770: {
771:   PetscReal sum = 0.0;
772:   PetscInt  d;
773:   for (d = 0; d < dim; ++d) sum += PetscRealPart(x[d]) * y[d];
774:   return sum;
775: }

777: static inline PetscReal DMPlex_DotRealD_Internal(PetscInt dim, const PetscReal *x, const PetscReal *y)
778: {
779:   PetscReal sum = 0.0;
780:   PetscInt  d;
781:   for (d = 0; d < dim; ++d) sum += x[d] * y[d];
782:   return sum;
783: }

785: static inline PetscReal DMPlex_NormD_Internal(PetscInt dim, const PetscReal *x)
786: {
787:   PetscReal sum = 0.0;
788:   PetscInt  d;
789:   for (d = 0; d < dim; ++d) sum += x[d] * x[d];
790:   return PetscSqrtReal(sum);
791: }

793: static inline PetscReal DMPlex_DistD_Internal(PetscInt dim, const PetscScalar *x, const PetscScalar *y)
794: {
795:   PetscReal sum = 0.0;
796:   PetscInt  d;
797:   for (d = 0; d < dim; ++d) sum += PetscRealPart(PetscConj(x[d] - y[d]) * (x[d] - y[d]));
798:   return PetscSqrtReal(sum);
799: }

801: PETSC_INTERN PetscErrorCode DMPlexGetPointDualSpaceFEM(DM, PetscInt, PetscInt, PetscDualSpace *);
802: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPoint_Internal(PetscSection, PetscBool, PetscInt, PetscInt, PetscInt *, PetscBool, const PetscInt[], const PetscInt[], PetscInt[]);
803: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPointFields_Internal(PetscSection, PetscBool, PetscInt, PetscInt, PetscInt[], PetscBool, const PetscInt ***, PetscInt, const PetscInt[], PetscInt[]);
804: PETSC_INTERN PetscErrorCode DMPlexGetTransitiveClosure_Internal(DM, PetscInt, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
805: PETSC_INTERN PetscErrorCode DMPlexVecGetOrientedClosure_Internal(DM, PetscSection, PetscBool, Vec, PetscInt, PetscInt, PetscInt *, PetscScalar *[]);
806: PETSC_INTERN PetscErrorCode DMPlexMatSetClosure_Internal(DM, PetscSection, PetscSection, PetscBool, Mat, PetscInt, const PetscScalar[], InsertMode);

808: PETSC_EXTERN PetscErrorCode DMPlexGetAllCells_Internal(DM, IS *);
809: PETSC_EXTERN PetscErrorCode DMPlexGetAllFaces_Internal(DM, IS *);
810: PETSC_EXTERN PetscErrorCode DMSNESGetFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
811: PETSC_EXTERN PetscErrorCode DMSNESRestoreFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
812: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Patch_Internal(DM, PetscSection, IS, PetscReal, Vec, Vec, Vec, void *);
813: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Patch_Internal(DM, PetscSection, PetscSection, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
814: PETSC_INTERN PetscErrorCode DMCreateSubDomainDM_Plex(DM, DMLabel, PetscInt, IS *, DM *);
815: PETSC_INTERN PetscErrorCode DMPlexBasisTransformPoint_Internal(DM, DM, Vec, PetscInt, PetscBool[], PetscBool, PetscScalar *);
816: PETSC_EXTERN PetscErrorCode DMPlexBasisTransformPointTensor_Internal(DM, DM, Vec, PetscInt, PetscBool, PetscInt, PetscScalar *);
817: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApplyReal_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscReal *, PetscReal *, void *);
818: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApply_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscScalar *, PetscScalar *, void *);
819: PETSC_INTERN PetscErrorCode DMCreateNeumannOverlap_Plex(DM, IS *, Mat *, PetscErrorCode (**)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void *), void **);
820: PETSC_INTERN PetscErrorCode DMPlexMarkBoundaryFaces_Internal(DM, PetscInt, PetscInt, DMLabel, PetscBool);
821: PETSC_INTERN PetscErrorCode DMPlexDistributeOverlap_Internal(DM, PetscInt, MPI_Comm, const char *, PetscSF *, DM *);

823: PETSC_INTERN PetscErrorCode DMPlexMarkSubmesh_Interpolated(DM, DMLabel, PetscInt, PetscBool, PetscBool, DMLabel, DM);

825: PETSC_INTERN PetscErrorCode DMPeriodicCoordinateSetUp_Internal(DM);

827: /* Functions in the vtable */
828: PETSC_INTERN PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
829: PETSC_INTERN PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat);
830: PETSC_INTERN PetscErrorCode DMCreateMassMatrix_Plex(DM dmCoarse, DM dmFine, Mat *mat);
831: PETSC_INTERN PetscErrorCode DMCreateMassMatrixLumped_Plex(DM, Vec *, Vec *);
832: PETSC_INTERN PetscErrorCode DMCreateLocalSection_Plex(DM dm);
833: PETSC_INTERN PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm);
834: PETSC_INTERN PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J);
835: PETSC_INTERN PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
836: PETSC_INTERN PetscErrorCode DMCreateCoordinateField_Plex(DM dm, DMField *field);
837: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
838: PETSC_INTERN PetscErrorCode DMSetUp_Plex(DM dm);
839: PETSC_INTERN PetscErrorCode DMDestroy_Plex(DM dm);
840: PETSC_INTERN PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
841: PETSC_INTERN PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
842: PETSC_INTERN PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm);
843: PETSC_INTERN PetscErrorCode DMCreateSuperDM_Plex(DM dms[], PetscInt len, IS **is, DM *superdm);
844: PETSC_INTERN PetscErrorCode DMCreateDomainDecompositionScatters_Plex(DM, PetscInt, DM *, VecScatter **, VecScatter **, VecScatter **);
845: PETSC_INTERN PetscErrorCode DMCreateDomainDecomposition_Plex(DM, PetscInt *, char ***, IS **, IS **, DM **);
846: PETSC_INTERN PetscErrorCode DMCreateSectionPermutation_Plex(DM dm, IS *permutation, PetscBT *blockStarts);

848: // Coordinate mapping functions
849: PETSC_INTERN void coordMap_identity(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
850: PETSC_INTERN void coordMap_shear(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
851: PETSC_INTERN void coordMap_flare(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
852: PETSC_INTERN void coordMap_annulus(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
853: PETSC_INTERN void coordMap_shell(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);

855: PETSC_EXTERN PetscErrorCode DMSnapToGeomModel_EGADS(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
856: PETSC_EXTERN PetscErrorCode DMSnapToGeomModel_EGADSLite(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);