Actual source code: dmpleximpl.h

petsc-master 2020-01-21
Report Typos and Errors
  1: #if !defined(_PLEXIMPL_H)
  2: #define _PLEXIMPL_H

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

 10: PETSC_EXTERN PetscLogEvent DMPLEX_Interpolate;
 11: PETSC_EXTERN PetscLogEvent DMPLEX_Partition;
 12: PETSC_EXTERN PetscLogEvent DMPLEX_PartSelf;
 13: PETSC_EXTERN PetscLogEvent DMPLEX_PartLabelInvert;
 14: PETSC_EXTERN PetscLogEvent DMPLEX_PartLabelCreateSF;
 15: PETSC_EXTERN PetscLogEvent DMPLEX_PartStratSF;
 16: PETSC_EXTERN PetscLogEvent DMPLEX_CreatePointSF;
 17: PETSC_EXTERN PetscLogEvent DMPLEX_Distribute;
 18: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeCones;
 19: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeLabels;
 20: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeSF;
 21: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeOverlap;
 22: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeField;
 23: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeData;
 24: PETSC_EXTERN PetscLogEvent DMPLEX_Migrate;
 25: PETSC_EXTERN PetscLogEvent DMPLEX_InterpolateSF;
 26: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalToNaturalBegin;
 27: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalToNaturalEnd;
 28: PETSC_EXTERN PetscLogEvent DMPLEX_NaturalToGlobalBegin;
 29: PETSC_EXTERN PetscLogEvent DMPLEX_NaturalToGlobalEnd;
 30: PETSC_EXTERN PetscLogEvent DMPLEX_Stratify;
 31: PETSC_EXTERN PetscLogEvent DMPLEX_Symmetrize;
 32: PETSC_EXTERN PetscLogEvent DMPLEX_Preallocate;
 33: PETSC_EXTERN PetscLogEvent DMPLEX_ResidualFEM;
 34: PETSC_EXTERN PetscLogEvent DMPLEX_JacobianFEM;
 35: PETSC_EXTERN PetscLogEvent DMPLEX_InterpolatorFEM;
 36: PETSC_EXTERN PetscLogEvent DMPLEX_InjectorFEM;
 37: PETSC_EXTERN PetscLogEvent DMPLEX_IntegralFEM;
 38: PETSC_EXTERN PetscLogEvent DMPLEX_CreateGmsh;
 39: PETSC_EXTERN PetscLogEvent DMPLEX_RebalanceSharedPoints;
 40: PETSC_EXTERN PetscLogEvent DMPLEX_CreateFromFile;
 41: PETSC_EXTERN PetscLogEvent DMPLEX_CreateFromCellList;
 42: PETSC_EXTERN PetscLogEvent DMPLEX_CreateFromCellList_Coordinates;

 44: PETSC_EXTERN PetscBool      PetscPartitionerRegisterAllCalled;
 45: PETSC_EXTERN PetscErrorCode PetscPartitionerRegisterAll(void);

 47: PETSC_EXTERN const char * const CellRefiners[];
 48: typedef enum {REFINER_NOOP = 0,
 49:               REFINER_SIMPLEX_1D,
 50:               REFINER_SIMPLEX_2D,
 51:               REFINER_HYBRID_SIMPLEX_2D,
 52:               REFINER_SIMPLEX_TO_HEX_2D,
 53:               REFINER_HYBRID_SIMPLEX_TO_HEX_2D,
 54:               REFINER_HEX_2D,
 55:               REFINER_HYBRID_HEX_2D,
 56:               REFINER_SIMPLEX_3D,
 57:               REFINER_HYBRID_SIMPLEX_3D,
 58:               REFINER_SIMPLEX_TO_HEX_3D,
 59:               REFINER_HYBRID_SIMPLEX_TO_HEX_3D,
 60:               REFINER_HEX_3D,
 61:               REFINER_HYBRID_HEX_3D} CellRefiner;

 63: typedef struct _PetscPartitionerOps *PetscPartitionerOps;
 64: struct _PetscPartitionerOps {
 65:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,PetscPartitioner);
 66:   PetscErrorCode (*setup)(PetscPartitioner);
 67:   PetscErrorCode (*view)(PetscPartitioner,PetscViewer);
 68:   PetscErrorCode (*destroy)(PetscPartitioner);
 69:   PetscErrorCode (*partition)(PetscPartitioner, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscSection, PetscSection, PetscSection, IS *);
 70: };

 72: struct _p_PetscPartitioner {
 73:   PETSCHEADER(struct _PetscPartitionerOps);
 74:   void        *data;            /* Implementation object */
 75:   PetscInt    height;           /* Height of points to partition into non-overlapping subsets */
 76:   PetscInt    edgeCut;          /* The number of edge cut by the partition */
 77:   PetscReal   balance;          /* The maximum partition size divided by the minimum size */
 78:   PetscViewer viewer;
 79:   PetscViewer viewerGraph;
 80:   PetscBool   viewGraph;
 81:   PetscBool   noGraph;          /* if true, the partitioner does not need the connectivity graph, only the number of local vertices */
 82:   PetscBool   usevwgt;          /* if true, the partitioner looks at the local section vertSection to weight the vertices of the graph */
 83: };

 85: typedef struct {
 86:   PetscInt dummy;
 87: } PetscPartitioner_Chaco;

 89: typedef struct {
 90:   MPI_Comm  pcomm;
 91:   PetscInt  ptype;
 92:   PetscReal imbalanceRatio;
 93:   PetscInt  debugFlag;
 94:   PetscInt  randomSeed;
 95: } PetscPartitioner_ParMetis;

 97: typedef struct {
 98:   MPI_Comm  pcomm;
 99:   PetscInt  strategy;
100:   PetscReal imbalance;
101: } PetscPartitioner_PTScotch;

103: static const char *const
104: PTScotchStrategyList[] = {
105:   "DEFAULT",
106:   "QUALITY",
107:   "SPEED",
108:   "BALANCE",
109:   "SAFETY",
110:   "SCALABILITY",
111:   "RECURSIVE",
112:   "REMAP"
113: };

115: typedef struct {
116:   PetscSection section;   /* Sizes for each partition */
117:   IS           partition; /* Points in each partition */
118:   PetscBool    random;    /* Flag for a random partition */
119: } PetscPartitioner_Shell;

121: typedef struct {
122:   PetscInt dummy;
123: } PetscPartitioner_Simple;

125: typedef struct {
126:   PetscInt dummy;
127: } PetscPartitioner_Gather;

129: /* Utility struct to store the contents of a Fluent file in memory */
130: typedef struct {
131:   int          index;    /* Type of section */
132:   unsigned int zoneID;
133:   unsigned int first;
134:   unsigned int last;
135:   int          type;
136:   int          nd;       /* Either ND or element-type */
137:   void        *data;
138: } FluentSection;

140: struct _PetscGridHash {
141:   PetscInt     dim;
142:   PetscReal    lower[3];    /* The lower-left corner */
143:   PetscReal    upper[3];    /* The upper-right corner */
144:   PetscReal    extent[3];   /* The box size */
145:   PetscReal    h[3];        /* The subbox size */
146:   PetscInt     n[3];        /* The number of subboxes */
147:   PetscSection cellSection; /* Offsets for cells in each subbox*/
148:   IS           cells;       /* List of cells in each subbox */
149:   DMLabel      cellsSparse; /* Sparse storage for cell map */
150: };

152: /* Point Numbering in Plex:

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

156:    First Stratum:  Cells [height 0]
157:    Second Stratum: Vertices [depth 0]
158:    Third Stratum:  Faces [height 1]
159:    Fourth Stratum: Edges [depth 1]

161:    We do this so that the numbering of a cell-vertex mesh does not change after interpolation. Within a given stratum,
162:    we allow additional segregation of points. When hybrid, or prismatic, points are present, the first hybrid point in
163:    a stratum is indicated by hybridPointMax[depth]. In addition, we allow ghost cells to be defined for use in finite
164:    volume methods, and these do not have full cones. These cells occur after any hybrid cells, and this division is
165:    indicated by ghostCellStart.
166: */
167: typedef struct {
168:   PetscInt             refct;

170:   PetscSection         coneSection;       /* Layout of cones (inedges for DAG) */
171:   PetscInt             maxConeSize;       /* Cached for fast lookup */
172:   PetscInt            *cones;             /* Cone for each point */
173:   PetscInt            *coneOrientations;  /* Orientation of each cone point, means cone traveral should start on point 'o', and if negative start on -(o+1) and go in reverse */
174:   PetscSection         supportSection;    /* Layout of cones (inedges for DAG) */
175:   PetscInt             maxSupportSize;    /* Cached for fast lookup */
176:   PetscInt            *supports;          /* Cone for each point */
177:   PetscBool            refinementUniform; /* Flag for uniform cell refinement */
178:   PetscReal            refinementLimit;   /* Maximum volume for refined cell */
179:   PetscErrorCode     (*refinementFunc)(const PetscReal [], PetscReal *); /* Function giving the maximum volume for refined cell */
180:   PetscInt             hybridPointMax[8]; /* Allow segregation of some points, each dimension has a divider (used in VTK output and refinement) */
181:   PetscInt             overlap;           /* Overlap of the partitions as passed to DMPlexDistribute() or DMPlexDistributeOverlap() */
182:   PetscInt             ghostCellStart;    /* The first ghost cell (for FV BC) or -1 */
183:   DMPlexInterpolatedFlag interpolated;
184:   DMPlexInterpolatedFlag interpolatedCollective;

186:   PetscInt            *facesTmp;          /* Work space for faces operation */

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

191:   /* Generation */
192:   char                *tetgenOpts;
193:   char                *triangleOpts;
194:   PetscPartitioner     partitioner;
195:   PetscBool            partitionBalance;  /* Evenly divide partition overlap when distributing */
196:   PetscBool            remeshBd;

198:   /* Submesh */
199:   DMLabel              subpointMap;       /* Label each original mesh point in the submesh with its depth, subpoint are the implicit numbering */

201:   /* Labels and numbering */
202:   PetscObjectState     depthState;        /* State of depth label, so that we can determine if a user changes it */
203:   PetscObjectState     celltypeState;     /* State of celltype label, so that we can determine if a user changes it */
204:   IS                   globalVertexNumbers;
205:   IS                   globalCellNumbers;

207:   /* Constraints */
208:   PetscSection         anchorSection;      /* maps constrained points to anchor points */
209:   IS                   anchorIS;           /* anchors indexed by the above section */
210:   PetscErrorCode     (*createanchors)(DM); /* automatically compute anchors (probably from tree constraints) */
211:   PetscErrorCode     (*computeanchormatrix)(DM,PetscSection,PetscSection,Mat);

213:   /* Tree: automatically construct constraints for hierarchically non-conforming meshes */
214:   PetscSection         parentSection;     /* dof == 1 if point has parent */
215:   PetscInt            *parents;           /* point to parent */
216:   PetscInt            *childIDs;          /* point to child ID */
217:   PetscSection         childSection;      /* inverse of parent section */
218:   PetscInt            *children;          /* point to children */
219:   DM                   referenceTree;     /* reference tree to which child ID's refer */
220:   PetscErrorCode      (*getchildsymmetry)(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*);

222:   /* MATIS support */
223:   PetscSection         subdomainSection;

225:   /* Adjacency */
226:   PetscBool            useAnchors;        /* Replace constrained points with their anchors in adjacency lists */
227:   PetscErrorCode      (*useradjacency)(DM,PetscInt,PetscInt*,PetscInt[],void*); /* User callback for adjacency */
228:   void                *useradjacencyctx;  /* User context for callback */

230:   /* Projection */
231:   PetscInt             maxProjectionHeight; /* maximum height of cells used in DMPlexProject functions */

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 */

242:   /* Debugging */
243:   PetscBool            printSetValues;
244:   PetscInt             printFEM;
245:   PetscInt             printL2;
246:   PetscReal            printTol;
247: } DM_Plex;

249: PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll_VTU(DM,PetscViewer);
250: PETSC_EXTERN PetscErrorCode VecView_Plex_Local(Vec,PetscViewer);
251: PETSC_EXTERN PetscErrorCode VecView_Plex_Native(Vec,PetscViewer);
252: PETSC_EXTERN PetscErrorCode VecView_Plex(Vec,PetscViewer);
253: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Local(Vec,PetscViewer);
254: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Native(Vec,PetscViewer);
255: PETSC_EXTERN PetscErrorCode VecLoad_Plex(Vec,PetscViewer);
256: PETSC_INTERN PetscErrorCode DMPlexGetFieldType_Internal(DM, PetscSection, PetscInt, PetscInt *, PetscInt *, PetscViewerVTKFieldType *);
257: PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM,PetscViewer);
258: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject,PetscViewer);
259: #if defined(PETSC_HAVE_HDF5)
260: PETSC_EXTERN PetscErrorCode VecView_Plex_Local_HDF5(Vec, PetscViewer);
261: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5(Vec, PetscViewer);
262: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5(Vec, PetscViewer);
263: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5_Native(Vec, PetscViewer);
264: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5_Native(Vec, PetscViewer);
265: PETSC_EXTERN PetscErrorCode DMPlexView_HDF5(DM, PetscViewer);
266: PETSC_EXTERN PetscErrorCode DMPlexLoad_HDF5(DM, PetscViewer);
267: #endif

269: PETSC_INTERN PetscErrorCode DMPlexClosurePoints_Private(DM,PetscInt,const PetscInt[],IS*);
270: PETSC_INTERN PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *, DM);
271: PETSC_INTERN PetscErrorCode DMCoarsen_Plex(DM, MPI_Comm, DM *);
272: PETSC_INTERN PetscErrorCode DMCoarsenHierarchy_Plex(DM, PetscInt, DM []);
273: PETSC_INTERN PetscErrorCode DMRefine_Plex(DM, MPI_Comm, DM *);
274: PETSC_INTERN PetscErrorCode DMRefineHierarchy_Plex(DM, PetscInt, DM []);
275: PETSC_INTERN PetscErrorCode DMAdaptLabel_Plex(DM, DMLabel, DM *);
276: PETSC_INTERN PetscErrorCode DMAdaptMetric_Plex(DM, Vec, DMLabel, DM *);
277: PETSC_INTERN PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
278: PETSC_INTERN PetscErrorCode DMProjectFunctionLocal_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
279: 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);
280: 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);
281: 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);
282: PETSC_INTERN PetscErrorCode DMComputeL2Diff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
283: PETSC_INTERN PetscErrorCode DMComputeL2GradientDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[], const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,const PetscReal [],PetscReal *);
284: PETSC_INTERN PetscErrorCode DMComputeL2FieldDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
285: PETSC_INTERN PetscErrorCode DMLocatePoints_Plex(DM, Vec, DMPointLocationType, PetscSF);

287: PETSC_INTERN PetscErrorCode DMPlexBuildFromCellList_Internal(DM, PetscInt, PetscInt, PetscInt, PetscInt, const int[], PetscBool);
288: PETSC_INTERN PetscErrorCode DMPlexBuildFromCellList_Parallel_Internal(DM, PetscInt, PetscInt, PetscInt, PetscInt, const int[], PetscBool, PetscSF *);
289: PETSC_INTERN PetscErrorCode DMPlexBuildCoordinates_Internal(DM, PetscInt, PetscInt, PetscInt, const double[]);
290: PETSC_INTERN PetscErrorCode DMPlexBuildCoordinates_Parallel_Internal(DM, PetscInt, PetscInt, PetscInt, PetscSF, const PetscReal[]);
291: PETSC_INTERN PetscErrorCode DMPlexLoadLabels_HDF5_Internal(DM, PetscViewer);
292: PETSC_INTERN PetscErrorCode DMPlexView_HDF5_Internal(DM, PetscViewer);
293: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Internal(DM, PetscViewer);
294: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Xdmf_Internal(DM, PetscViewer);
295: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Internal(Vec, PetscViewer);
296: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Native_Internal(Vec, PetscViewer);
297: PETSC_INTERN PetscErrorCode VecView_Plex_Local_HDF5_Internal(Vec, PetscViewer);
298: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Internal(Vec, PetscViewer);
299: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Native_Internal(Vec, PetscViewer);
300: /* TODO Make these INTERN */
301: PETSC_EXTERN PetscErrorCode DMPlexView_ExodusII_Internal(DM, int, PetscInt);
302: PETSC_EXTERN PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec, int, int);
303: PETSC_EXTERN PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec, int, int);
304: PETSC_EXTERN PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec, int, int);
305: PETSC_EXTERN PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec, int, int);
306: PETSC_INTERN PetscErrorCode DMPlexVTKGetCellType_Internal(DM,PetscInt,PetscInt,PetscInt*);
307: PETSC_INTERN PetscErrorCode DMPlexGetAdjacency_Internal(DM,PetscInt,PetscBool,PetscBool,PetscBool,PetscInt*,PetscInt*[]);
308: PETSC_INTERN PetscErrorCode DMPlexGetFaces_Internal(DM,PetscInt,PetscInt*,PetscInt*,const PetscInt*[]);
309: PETSC_INTERN PetscErrorCode DMPlexGetRawFaces_Internal(DM,DMPolytopeType,const PetscInt[], PetscInt*,PetscInt*,const PetscInt*[]);
310: PETSC_INTERN PetscErrorCode DMPlexRestoreFaces_Internal(DM,PetscInt,PetscInt*,PetscInt*,const PetscInt*[]);
311: PETSC_INTERN PetscErrorCode DMPlexRefineUniform_Internal(DM,CellRefiner,DM*);
312: PETSC_INTERN PetscErrorCode DMPlexGetCellRefiner_Internal(DM,CellRefiner*);
313: PETSC_INTERN PetscErrorCode CellRefinerGetAffineTransforms_Internal(CellRefiner, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[]);
314: PETSC_INTERN PetscErrorCode CellRefinerGetAffineFaceTransforms_Internal(CellRefiner, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[], PetscReal *[]);
315: PETSC_INTERN PetscErrorCode CellRefinerRestoreAffineTransforms_Internal(CellRefiner, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[]);
316: PETSC_INTERN PetscErrorCode CellRefinerInCellTest_Internal(CellRefiner, const PetscReal[], PetscBool *);
317: PETSC_INTERN PetscErrorCode DMPlexInvertCell_Internal(PetscInt, PetscInt, PetscInt[]);
318: PETSC_INTERN PetscErrorCode DMPlexVecSetFieldClosure_Internal(DM, PetscSection, Vec, PetscBool[], PetscInt, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
319: PETSC_INTERN PetscErrorCode DMPlexProjectConstraints_Internal(DM, Vec, Vec);
320: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_SetTree(DM, PetscSection, PetscInt[], PetscInt[]);
321: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_Union(DM,DM,const char *,DM*);
322: PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorTree(DM,DM,PetscSF,PetscInt *,Mat);
323: PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorTree(DM,DM,PetscSF,PetscInt *,Mat);
324: PETSC_EXTERN PetscErrorCode DMPlexAnchorsModifyMat(DM,PetscSection,PetscInt,PetscInt,const PetscInt[],const PetscInt ***,const PetscScalar[],PetscInt*,PetscInt*,PetscInt*[],PetscScalar*[],PetscInt[],PetscBool);
325: PETSC_EXTERN PetscErrorCode indicesPoint_private(PetscSection,PetscInt,PetscInt,PetscInt *,PetscBool,PetscInt,PetscInt []);
326: PETSC_EXTERN PetscErrorCode indicesPointFields_private(PetscSection,PetscInt,PetscInt,PetscInt [],PetscBool,PetscInt,PetscInt []);
327: PETSC_INTERN PetscErrorCode DMPlexLocatePoint_Internal(DM,PetscInt,const PetscScalar [],PetscInt,PetscInt *);
328: /* these two are PETSC_EXTERN just because of src/dm/impls/plex/examples/tests/ex18.c */
329: PETSC_EXTERN PetscErrorCode DMPlexOrientCell_Internal(DM,PetscInt,PetscInt,PetscBool);
330: PETSC_EXTERN PetscErrorCode DMPlexOrientInterface_Internal(DM);

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

335: PETSC_INTERN PetscErrorCode DMPlexCreateCellNumbering_Internal(DM, PetscBool, IS *);
336: PETSC_INTERN PetscErrorCode DMPlexCreateVertexNumbering_Internal(DM, PetscBool, IS *);
337: PETSC_INTERN PetscErrorCode DMPlexRefine_Internal(DM, DMLabel, DM *);
338: PETSC_INTERN PetscErrorCode DMPlexCoarsen_Internal(DM, DMLabel, DM *);
339: PETSC_INTERN PetscErrorCode DMCreateMatrix_Plex(DM, Mat*);

341: PETSC_INTERN PetscErrorCode DMPlexGetOverlap_Plex(DM, PetscInt *);

343: /* invert dihedral symmetry: return a^-1,
344:  * using the representation described in
345:  * DMPlexGetConeOrientation() */
346: PETSC_STATIC_INLINE PetscInt DihedralInvert(PetscInt N, PetscInt a)
347: {
348:   return (a <= 0) ? a : (N - a);
349: }

351: /* invert dihedral symmetry: return b * a,
352:  * using the representation described in
353:  * DMPlexGetConeOrientation() */
354: PETSC_STATIC_INLINE PetscInt DihedralCompose(PetscInt N, PetscInt a, PetscInt b)
355: {
356:   if (!N) return 0;
357:   return  (a >= 0) ?
358:          ((b >= 0) ? ((a + b) % N) : -(((a - b - 1) % N) + 1)) :
359:          ((b >= 0) ? -(((N - b - a - 1) % N) + 1) : ((N + b - a) % N));
360: }

362: /* swap dihedral symmetries: return b * a^-1,
363:  * using the representation described in
364:  * DMPlexGetConeOrientation() */
365: PETSC_STATIC_INLINE PetscInt DihedralSwap(PetscInt N, PetscInt a, PetscInt b)
366: {
367:   return DihedralCompose(N,DihedralInvert(N,a),b);
368: }

370: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Internal(DM, IS , PetscReal, Vec, Vec, PetscReal, Vec, void *);
371: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Internal(DM, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
372: PETSC_EXTERN PetscErrorCode DMPlexReconstructGradients_Internal(DM, PetscFV, PetscInt, PetscInt, Vec, Vec, Vec, Vec);

374: /* Matvec with A in row-major storage, x and y can be aliased */
375: PETSC_STATIC_INLINE void DMPlex_Mult2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
376: {
377:   PetscScalar z[2];
378:   z[0] = x[0]; z[1] = x[ldx];
379:   y[0]   = A[0]*z[0] + A[1]*z[1];
380:   y[ldx] = A[2]*z[0] + A[3]*z[1];
381:   (void)PetscLogFlops(6.0);
382: }
383: PETSC_STATIC_INLINE void DMPlex_Mult3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
384: {
385:   PetscScalar z[3];
386:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
387:   y[0]     = A[0]*z[0] + A[1]*z[1] + A[2]*z[2];
388:   y[ldx]   = A[3]*z[0] + A[4]*z[1] + A[5]*z[2];
389:   y[ldx*2] = A[6]*z[0] + A[7]*z[1] + A[8]*z[2];
390:   (void)PetscLogFlops(15.0);
391: }
392: PETSC_STATIC_INLINE void DMPlex_MultTranspose2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
393: {
394:   PetscScalar z[2];
395:   z[0] = x[0]; z[1] = x[ldx];
396:   y[0]   = A[0]*z[0] + A[2]*z[1];
397:   y[ldx] = A[1]*z[0] + A[3]*z[1];
398:   (void)PetscLogFlops(6.0);
399: }
400: PETSC_STATIC_INLINE void DMPlex_MultTranspose3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
401: {
402:   PetscScalar z[3];
403:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
404:   y[0]     = A[0]*z[0] + A[3]*z[1] + A[6]*z[2];
405:   y[ldx]   = A[1]*z[0] + A[4]*z[1] + A[7]*z[2];
406:   y[ldx*2] = A[2]*z[0] + A[5]*z[1] + A[8]*z[2];
407:   (void)PetscLogFlops(15.0);
408: }
409: PETSC_STATIC_INLINE void DMPlex_Mult2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
410: {
411:   PetscScalar z[2];
412:   z[0] = x[0]; z[1] = x[ldx];
413:   y[0]   = A[0]*z[0] + A[1]*z[1];
414:   y[ldx] = A[2]*z[0] + A[3]*z[1];
415:   (void)PetscLogFlops(6.0);
416: }
417: PETSC_STATIC_INLINE void DMPlex_Mult3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
418: {
419:   PetscScalar z[3];
420:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
421:   y[0]     = A[0]*z[0] + A[1]*z[1] + A[2]*z[2];
422:   y[ldx]   = A[3]*z[0] + A[4]*z[1] + A[5]*z[2];
423:   y[ldx*2] = A[6]*z[0] + A[7]*z[1] + A[8]*z[2];
424:   (void)PetscLogFlops(15.0);
425: }
426: PETSC_STATIC_INLINE void DMPlex_MultTranspose2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
427: {
428:   PetscScalar z[2];
429:   z[0] = x[0]; z[1] = x[ldx];
430:   y[0]   = A[0]*z[0] + A[2]*z[1];
431:   y[ldx] = A[1]*z[0] + A[3]*z[1];
432:   (void)PetscLogFlops(6.0);
433: }
434: PETSC_STATIC_INLINE void DMPlex_MultTranspose3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
435: {
436:   PetscScalar z[3];
437:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
438:   y[0]     = A[0]*z[0] + A[3]*z[1] + A[6]*z[2];
439:   y[ldx]   = A[1]*z[0] + A[4]*z[1] + A[7]*z[2];
440:   y[ldx*2] = A[2]*z[0] + A[5]*z[1] + A[8]*z[2];
441:   (void)PetscLogFlops(15.0);
442: }

444: PETSC_STATIC_INLINE void DMPlex_MatMult2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
445: {
446:   PetscInt j;
447:   for (j = 0; j < n; ++j) {
448:     PetscScalar z[2];
449:     z[0] = B[0+j]; z[1] = B[1*ldb+j];
450:     DMPlex_Mult2D_Internal(A, 1, z, z);
451:     C[0+j] = z[0]; C[1*ldb+j] = z[1];
452:   }
453:   (void)PetscLogFlops(8.0*n);
454: }
455: PETSC_STATIC_INLINE void DMPlex_MatMult3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
456: {
457:   PetscInt j;
458:   for (j = 0; j < n; ++j) {
459:     PetscScalar z[3];
460:     z[0] = B[0+j]; z[1] = B[1*ldb+j]; z[2] = B[2*ldb+j];
461:     DMPlex_Mult3D_Internal(A, 1, z, z);
462:     C[0+j] = z[0]; C[1*ldb+j] = z[1]; C[2*ldb+j] = z[2];
463:   }
464:   (void)PetscLogFlops(8.0*n);
465: }
466: PETSC_STATIC_INLINE void DMPlex_MatMultTranspose2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
467: {
468:   PetscInt j;
469:   for (j = 0; j < n; ++j) {
470:     PetscScalar z[2];
471:     z[0] = B[0+j]; z[1] = B[1*ldb+j];
472:     DMPlex_MultTranspose2D_Internal(A, 1, z, z);
473:     C[0+j] = z[0]; C[1*ldb+j] = z[1];
474:   }
475:   (void)PetscLogFlops(8.0*n);
476: }
477: PETSC_STATIC_INLINE void DMPlex_MatMultTranspose3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
478: {
479:   PetscInt j;
480:   for (j = 0; j < n; ++j) {
481:     PetscScalar z[3];
482:     z[0] = B[0+j]; z[1] = B[1*ldb+j]; z[2] = B[2*ldb+j];
483:     DMPlex_MultTranspose3D_Internal(A, 1, z, z);
484:     C[0+j] = z[0]; C[1*ldb+j] = z[1]; C[2*ldb+j] = z[2];
485:   }
486:   (void)PetscLogFlops(8.0*n);
487: }

489: PETSC_STATIC_INLINE void DMPlex_MatMultLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
490: {
491:   PetscInt j;
492:   for (j = 0; j < m; ++j) {
493:     DMPlex_MultTranspose2D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
494:   }
495:   (void)PetscLogFlops(8.0*m);
496: }
497: PETSC_STATIC_INLINE void DMPlex_MatMultLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
498: {
499:   PetscInt j;
500:   for (j = 0; j < m; ++j) {
501:     DMPlex_MultTranspose3D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
502:   }
503:   (void)PetscLogFlops(8.0*m);
504: }
505: PETSC_STATIC_INLINE void DMPlex_MatMultTransposeLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
506: {
507:   PetscInt j;
508:   for (j = 0; j < m; ++j) {
509:     DMPlex_Mult2D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
510:   }
511:   (void)PetscLogFlops(8.0*m);
512: }
513: PETSC_STATIC_INLINE void DMPlex_MatMultTransposeLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
514: {
515:   PetscInt j;
516:   for (j = 0; j < m; ++j) {
517:     DMPlex_Mult3D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
518:   }
519:   (void)PetscLogFlops(8.0*m);
520: }

522: PETSC_STATIC_INLINE void DMPlex_Transpose2D_Internal(PetscScalar A[])
523: {
524:   PetscScalar tmp;
525:   tmp = A[1]; A[1] = A[2]; A[2] = tmp;
526: }
527: PETSC_STATIC_INLINE void DMPlex_Transpose3D_Internal(PetscScalar A[])
528: {
529:   PetscScalar tmp;
530:   tmp = A[1]; A[1] = A[3]; A[3] = tmp;
531:   tmp = A[2]; A[2] = A[6]; A[6] = tmp;
532:   tmp = A[5]; A[5] = A[7]; A[7] = tmp;
533: }

535: PETSC_STATIC_INLINE void DMPlex_Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
536: {
537:   const PetscReal invDet = 1.0/detJ;

539:   invJ[0] =  invDet*J[3];
540:   invJ[1] = -invDet*J[1];
541:   invJ[2] = -invDet*J[2];
542:   invJ[3] =  invDet*J[0];
543:   (void)PetscLogFlops(5.0);
544: }

546: PETSC_STATIC_INLINE void DMPlex_Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
547: {
548:   const PetscReal invDet = 1.0/detJ;

550:   invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
551:   invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
552:   invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
553:   invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
554:   invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
555:   invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
556:   invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
557:   invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
558:   invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
559:   (void)PetscLogFlops(37.0);
560: }

562: PETSC_STATIC_INLINE void DMPlex_Det2D_Internal(PetscReal *detJ, const PetscReal J[])
563: {
564:   *detJ = J[0]*J[3] - J[1]*J[2];
565:   (void)PetscLogFlops(3.0);
566: }

568: PETSC_STATIC_INLINE void DMPlex_Det3D_Internal(PetscReal *detJ, const PetscReal J[])
569: {
570:   *detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
571:            J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
572:            J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
573:   (void)PetscLogFlops(12.0);
574: }

576: PETSC_STATIC_INLINE void DMPlex_Det2D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
577: {
578:   *detJ = PetscRealPart(J[0])*PetscRealPart(J[3]) - PetscRealPart(J[1])*PetscRealPart(J[2]);
579:   (void)PetscLogFlops(3.0);
580: }

582: PETSC_STATIC_INLINE void DMPlex_Det3D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
583: {
584:   *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])) +
585:            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])) +
586:            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])));
587:   (void)PetscLogFlops(12.0);
588: }

590: PETSC_STATIC_INLINE void DMPlex_WaxpyD_Internal(PetscInt dim, PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w) {PetscInt d; for (d = 0; d < dim; ++d) w[d] = a*x[d] + y[d];}

592: PETSC_STATIC_INLINE PetscReal DMPlex_DotD_Internal(PetscInt dim, const PetscScalar *x, const PetscReal *y) {PetscReal sum = 0.0; PetscInt d; for (d = 0; d < dim; ++d) sum += PetscRealPart(x[d])*y[d]; return sum;}

594: PETSC_STATIC_INLINE PetscReal DMPlex_DotRealD_Internal(PetscInt dim, const PetscReal *x, const PetscReal *y) {PetscReal sum = 0.0; PetscInt d; for (d = 0; d < dim; ++d) sum += x[d]*y[d]; return sum;}

596: PETSC_STATIC_INLINE PetscReal DMPlex_NormD_Internal(PetscInt dim, const PetscReal *x) {PetscReal sum = 0.0; PetscInt d; for (d = 0; d < dim; ++d) sum += x[d]*x[d]; return PetscSqrtReal(sum);}

598: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Translate_Private(PetscInt ornt, PetscInt *start, PetscBool *reverse)
599: {
601:   *reverse = (ornt < 0) ? PETSC_TRUE : PETSC_FALSE;
602:   *start = *reverse ? -(ornt+1) : ornt;
603:   return(0);
604: }

606: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Combine_Private(PetscInt coneSize, PetscInt origStart, PetscBool origReverse, PetscInt rotateStart, PetscBool rotateReverse, PetscInt *newStart, PetscBool *newReverse)
607: {
609:   *newReverse = (origReverse == rotateReverse) ? PETSC_FALSE : PETSC_TRUE;
610:   *newStart = rotateReverse ? (coneSize + rotateStart - origStart) : (coneSize + origStart - rotateStart);
611:   *newStart %= coneSize;
612:   return(0);
613: }

615: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_TranslateBack_Private(PetscInt coneSize, PetscInt start, PetscBool reverse, PetscInt *ornt)
616: {
618:   if (coneSize < 3) {
619:     /* edges just get flipped if start == 1 regardless direction */
620:     *ornt = start ? -2 : 0;
621:   } else {
622:     *ornt = reverse ? -(start+1) : start;
623:   }
624:   return(0);
625: }

627: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Permute_Private(PetscInt n, const PetscInt arr[], PetscInt start, PetscBool reverse, PetscInt newarr[])
628: {
629:   PetscInt i;

632:   if (reverse) {for (i=0; i<n; i++) newarr[i] = arr[(n+start-i)%n];}
633:   else         {for (i=0; i<n; i++) newarr[i] = arr[(start+i)%n];}
634:   return(0);
635: }

637: PETSC_INTERN PetscErrorCode DMPlexGetPointDualSpaceFEM(DM,PetscInt,PetscInt,PetscDualSpace *);
638: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPoint_Internal(PetscSection,PetscInt,PetscInt,PetscInt *,PetscBool,const PetscInt[],const PetscInt[],PetscInt[]);
639: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPointFields_Internal(PetscSection,PetscInt,PetscInt,PetscInt[],PetscBool,const PetscInt***,PetscInt,const PetscInt[],PetscInt[]);
640: PETSC_INTERN PetscErrorCode DMPlexGetCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);
641: PETSC_INTERN PetscErrorCode DMPlexRestoreCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);

643: PETSC_EXTERN PetscErrorCode DMSNESGetFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
644: PETSC_EXTERN PetscErrorCode DMSNESRestoreFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
645: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Patch_Internal(DM, PetscSection, IS, PetscReal, Vec, Vec, Vec, void *);
646: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Patch_Internal(DM, PetscSection, PetscSection, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
647: PETSC_INTERN PetscErrorCode DMCreateSubDomainDM_Plex(DM,DMLabel,PetscInt,IS*,DM*);
648: PETSC_INTERN PetscErrorCode DMPlexBasisTransformPoint_Internal(DM, DM, Vec, PetscInt, PetscBool[], PetscBool, PetscScalar *);
649: PETSC_EXTERN PetscErrorCode DMPlexBasisTransformPointTensor_Internal(DM, DM, Vec, PetscInt, PetscBool, PetscInt, PetscScalar *);
650: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApplyReal_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscReal *, PetscReal *, void *);
651: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApply_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscScalar *, PetscScalar *, void *);
652: PETSC_INTERN PetscErrorCode DMCreateNeumannOverlap_Plex(DM, IS*, Mat*, PetscErrorCode (**)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*), void **);

654: #endif /* _PLEXIMPL_H */