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_CreateBoxSFC;
 47: PETSC_EXTERN PetscLogEvent DMPLEX_RebalanceSharedPoints;
 48: PETSC_EXTERN PetscLogEvent DMPLEX_CreateFromFile;
 49: PETSC_EXTERN PetscLogEvent DMPLEX_CreateFromOptions;
 50: PETSC_EXTERN PetscLogEvent DMPLEX_BuildFromCellList;
 51: PETSC_EXTERN PetscLogEvent DMPLEX_BuildCoordinatesFromCellList;
 52: PETSC_EXTERN PetscLogEvent DMPLEX_LocatePoints;
 53: PETSC_EXTERN PetscLogEvent DMPLEX_TopologyView;
 54: PETSC_EXTERN PetscLogEvent DMPLEX_DistributionView;
 55: PETSC_EXTERN PetscLogEvent DMPLEX_LabelsView;
 56: PETSC_EXTERN PetscLogEvent DMPLEX_CoordinatesView;
 57: PETSC_EXTERN PetscLogEvent DMPLEX_SectionView;
 58: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalVectorView;
 59: PETSC_EXTERN PetscLogEvent DMPLEX_LocalVectorView;
 60: PETSC_EXTERN PetscLogEvent DMPLEX_TopologyLoad;
 61: PETSC_EXTERN PetscLogEvent DMPLEX_DistributionLoad;
 62: PETSC_EXTERN PetscLogEvent DMPLEX_LabelsLoad;
 63: PETSC_EXTERN PetscLogEvent DMPLEX_CoordinatesLoad;
 64: PETSC_EXTERN PetscLogEvent DMPLEX_SectionLoad;
 65: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalVectorLoad;
 66: PETSC_EXTERN PetscLogEvent DMPLEX_LocalVectorLoad;
 67: PETSC_EXTERN PetscLogEvent DMPLEX_MetricEnforceSPD;
 68: PETSC_EXTERN PetscLogEvent DMPLEX_MetricNormalize;
 69: PETSC_EXTERN PetscLogEvent DMPLEX_MetricAverage;
 70: PETSC_EXTERN PetscLogEvent DMPLEX_MetricIntersection;
 71: PETSC_EXTERN PetscLogEvent DMPLEX_Generate;
 72: PETSC_EXTERN PetscLogEvent DMPLEX_Transform;
 73: PETSC_EXTERN PetscLogEvent DMPLEX_GetLocalOffsets;

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

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

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

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

123: /* Point Numbering in Plex:

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

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

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

138:   PetscSection coneSection;      /* Layout of cones (inedges for DAG) */
139:   PetscInt    *cones;            /* Cone for each point */
140:   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 */
141:   PetscSection supportSection;   /* Layout of cones (inedges for DAG) */
142:   PetscInt    *supports;         /* Cone for each point */

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

237:   /* Geometry */
238:   PetscReal     minradius;                        /* Minimum distance from cell centroid to face */
239:   PetscBool     useHashLocation;                  /* Use grid hashing for point location */
240:   PetscGridHash lbox;                             /* Local box for searching */
241:   void (*coordFunc)(PetscInt, PetscInt, PetscInt, /* Function used to remap newly introduced vertices */
242:                     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[]);

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

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

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

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

264: PETSC_INTERN PetscErrorCode DMPlexCopy_Internal(DM, PetscBool, PetscBool, DM);
265: PETSC_INTERN PetscErrorCode DMPlexReplace_Internal(DM, DM *);

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

312: PETSC_INTERN PetscErrorCode DMPlexVecGetClosure_Internal(DM, PetscSection, PetscBool, Vec, PetscInt, PetscInt *, PetscScalar *[]);
313: PETSC_INTERN PetscErrorCode DMPlexVecGetClosureAtDepth_Internal(DM, PetscSection, Vec, PetscInt, PetscInt, PetscInt *, PetscScalar *[]);
314: PETSC_INTERN PetscErrorCode DMPlexClosurePoints_Private(DM, PetscInt, const PetscInt[], IS *);
315: PETSC_INTERN PetscErrorCode DMSetFromOptions_NonRefinement_Plex(DM, PetscOptionItems *);
316: PETSC_INTERN PetscErrorCode DMSetFromOptions_Overlap_Plex(DM, PetscOptionItems *, PetscInt *);
317: PETSC_INTERN PetscErrorCode DMCoarsen_Plex(DM, MPI_Comm, DM *);
318: PETSC_INTERN PetscErrorCode DMCoarsenHierarchy_Plex(DM, PetscInt, DM[]);
319: PETSC_INTERN PetscErrorCode DMRefine_Plex(DM, MPI_Comm, DM *);
320: PETSC_INTERN PetscErrorCode DMRefineHierarchy_Plex(DM, PetscInt, DM[]);
321: PETSC_INTERN PetscErrorCode DMAdaptLabel_Plex(DM, Vec, DMLabel, DMLabel, DM *);
322: PETSC_INTERN PetscErrorCode DMExtrude_Plex(DM, PetscInt, DM *);
323: PETSC_INTERN PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
324: PETSC_INTERN PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
325: PETSC_INTERN PetscErrorCode DMProjectFunctionLocal_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
326: 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);
327: 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);
328: 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);
329: 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);
330: PETSC_INTERN PetscErrorCode DMComputeL2Diff_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
331: PETSC_INTERN PetscErrorCode DMComputeL2GradientDiff_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal[], PetscReal *);
332: PETSC_INTERN PetscErrorCode DMComputeL2FieldDiff_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
333: PETSC_INTERN PetscErrorCode DMLocatePoints_Plex(DM, Vec, DMPointLocationType, PetscSF);

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

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

366: PETSC_INTERN PetscErrorCode DMPlexInterpolateInPlace_Internal(DM);
367: PETSC_INTERN PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Internal(DM, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool);
368: PETSC_INTERN PetscErrorCode DMPlexMigrateIsoperiodicFaceSF_Internal(DM, DM, PetscSF);
369: PETSC_INTERN PetscErrorCode DMPlexCreateVertexNumbering_Internal(DM, PetscBool, IS *);
370: PETSC_INTERN PetscErrorCode DMPlexRefine_Internal(DM, Vec, DMLabel, DMLabel, DM *);
371: PETSC_INTERN PetscErrorCode DMPlexCoarsen_Internal(DM, Vec, DMLabel, DMLabel, DM *);
372: PETSC_INTERN PetscErrorCode DMCreateMatrix_Plex(DM, Mat *);

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

387: #if 1
388: static inline PetscInt DihedralInvert(PetscInt N, PetscInt a)
389: {
390:   return (a <= 0) ? a : (N - a);
391: }

393: static inline PetscInt DihedralCompose(PetscInt N, PetscInt a, PetscInt b)
394: {
395:   if (!N) return 0;
396:   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));
397: }

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

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

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

472: PETSC_EXTERN PetscErrorCode DMPlexComputeIntegral_Internal(DM, Vec, PetscInt, PetscInt, PetscScalar *, void *);
473: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Internal(DM, PetscFormKey, IS, PetscReal, Vec, Vec, PetscReal, Vec, void *);
474: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Hybrid_Internal(DM, PetscFormKey[], IS, PetscReal, Vec, Vec, PetscReal, Vec, void *);
475: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Internal(DM, PetscFormKey, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
476: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Hybrid_Internal(DM, PetscFormKey[], IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
477: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Action_Internal(DM, PetscFormKey, IS, PetscReal, PetscReal, Vec, Vec, Vec, Vec, void *);
478: PETSC_EXTERN PetscErrorCode DMPlexReconstructGradients_Internal(DM, PetscFV, PetscInt, PetscInt, Vec, Vec, Vec, Vec);

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

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

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

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

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

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

714: static inline void DMPlex_Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
715: {
716:   // Allow zero volume cells
717:   const PetscReal invDet = detJ == 0 ? 1.0 : (PetscReal)1.0 / detJ;

719:   invJ[0] = invDet * J[3];
720:   invJ[1] = -invDet * J[1];
721:   invJ[2] = -invDet * J[2];
722:   invJ[3] = invDet * J[0];
723:   (void)PetscLogFlops(5.0);
724: }

726: static inline void DMPlex_Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
727: {
728:   // Allow zero volume cells
729:   const PetscReal invDet = detJ == 0 ? 1.0 : (PetscReal)1.0 / detJ;

731:   invJ[0 * 3 + 0] = invDet * (J[1 * 3 + 1] * J[2 * 3 + 2] - J[1 * 3 + 2] * J[2 * 3 + 1]);
732:   invJ[0 * 3 + 1] = invDet * (J[0 * 3 + 2] * J[2 * 3 + 1] - J[0 * 3 + 1] * J[2 * 3 + 2]);
733:   invJ[0 * 3 + 2] = invDet * (J[0 * 3 + 1] * J[1 * 3 + 2] - J[0 * 3 + 2] * J[1 * 3 + 1]);
734:   invJ[1 * 3 + 0] = invDet * (J[1 * 3 + 2] * J[2 * 3 + 0] - J[1 * 3 + 0] * J[2 * 3 + 2]);
735:   invJ[1 * 3 + 1] = invDet * (J[0 * 3 + 0] * J[2 * 3 + 2] - J[0 * 3 + 2] * J[2 * 3 + 0]);
736:   invJ[1 * 3 + 2] = invDet * (J[0 * 3 + 2] * J[1 * 3 + 0] - J[0 * 3 + 0] * J[1 * 3 + 2]);
737:   invJ[2 * 3 + 0] = invDet * (J[1 * 3 + 0] * J[2 * 3 + 1] - J[1 * 3 + 1] * J[2 * 3 + 0]);
738:   invJ[2 * 3 + 1] = invDet * (J[0 * 3 + 1] * J[2 * 3 + 0] - J[0 * 3 + 0] * J[2 * 3 + 1]);
739:   invJ[2 * 3 + 2] = invDet * (J[0 * 3 + 0] * J[1 * 3 + 1] - J[0 * 3 + 1] * J[1 * 3 + 0]);
740:   (void)PetscLogFlops(37.0);
741: }

743: static inline void DMPlex_Det2D_Internal(PetscReal *detJ, const PetscReal J[])
744: {
745:   *detJ = J[0] * J[3] - J[1] * J[2];
746:   (void)PetscLogFlops(3.0);
747: }

749: static inline void DMPlex_Det3D_Internal(PetscReal *detJ, const PetscReal J[])
750: {
751:   *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]));
752:   (void)PetscLogFlops(12.0);
753: }

755: static inline void DMPlex_Det2D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
756: {
757:   *detJ = PetscRealPart(J[0]) * PetscRealPart(J[3]) - PetscRealPart(J[1]) * PetscRealPart(J[2]);
758:   (void)PetscLogFlops(3.0);
759: }

761: static inline void DMPlex_Det3D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
762: {
763:   *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])));
764:   (void)PetscLogFlops(12.0);
765: }

767: static inline void DMPlex_WaxpyD_Internal(PetscInt dim, PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w)
768: {
769:   PetscInt d;
770:   for (d = 0; d < dim; ++d) w[d] = a * x[d] + y[d];
771: }

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

781: static inline PetscReal DMPlex_DotRealD_Internal(PetscInt dim, const PetscReal *x, const PetscReal *y)
782: {
783:   PetscReal sum = 0.0;
784:   PetscInt  d;
785:   for (d = 0; d < dim; ++d) sum += x[d] * y[d];
786:   return sum;
787: }

789: static inline PetscReal DMPlex_NormD_Internal(PetscInt dim, const PetscReal *x)
790: {
791:   PetscReal sum = 0.0;
792:   PetscInt  d;
793:   for (d = 0; d < dim; ++d) sum += x[d] * x[d];
794:   return PetscSqrtReal(sum);
795: }

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

805: PETSC_INTERN PetscErrorCode DMPlexGetPointDualSpaceFEM(DM, PetscInt, PetscInt, PetscDualSpace *);
806: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPoint_Internal(PetscSection, PetscBool, PetscInt, PetscInt, PetscInt *, PetscBool, const PetscInt[], const PetscInt[], PetscInt[]);
807: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPointFields_Internal(PetscSection, PetscBool, PetscInt, PetscInt, PetscInt[], PetscBool, const PetscInt ***, PetscInt, const PetscInt[], PetscInt[]);
808: PETSC_INTERN PetscErrorCode DMPlexGetTransitiveClosure_Internal(DM, PetscInt, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
809: PETSC_INTERN PetscErrorCode DMPlexVecGetOrientedClosure_Internal(DM, PetscSection, PetscBool, Vec, PetscInt, PetscInt, PetscInt *, PetscScalar *[]);
810: PETSC_INTERN PetscErrorCode DMPlexMatSetClosure_Internal(DM, PetscSection, PetscSection, PetscBool, Mat, PetscInt, const PetscScalar[], InsertMode);

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

827: PETSC_INTERN PetscErrorCode DMPlexInterpolateFaces_Internal(DM, PetscInt, DM);

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

831: PETSC_INTERN PetscErrorCode DMPeriodicCoordinateSetUp_Internal(DM);

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

854: // Coordinate mapping functions
855: 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[]);
856: 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[]);
857: 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[]);
858: 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[]);
859: 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[]);

861: PETSC_EXTERN PetscErrorCode DMSnapToGeomModel_EGADS(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
862: PETSC_EXTERN PetscErrorCode DMSnapToGeomModel_EGADSLite(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);