Actual source code: dmpleximpl.h

petsc-master 2019-11-16
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;

 41: PETSC_EXTERN PetscBool      PetscPartitionerRegisterAllCalled;
 42: PETSC_EXTERN PetscErrorCode PetscPartitionerRegisterAll(void);

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

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

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

 81: typedef struct {
 82:   PetscInt dummy;
 83: } PetscPartitioner_Chaco;

 85: typedef struct {
 86:   PetscInt  ptype;
 87:   PetscReal imbalanceRatio;
 88:   PetscInt  debugFlag;
 89:   PetscInt  randomSeed;
 90: } PetscPartitioner_ParMetis;

 92: typedef struct {
 93:   PetscInt  strategy;
 94:   PetscReal imbalance;
 95: } PetscPartitioner_PTScotch;

 97: static const char *const
 98: PTScotchStrategyList[] = {
 99:   "DEFAULT",
100:   "QUALITY",
101:   "SPEED",
102:   "BALANCE",
103:   "SAFETY",
104:   "SCALABILITY",
105:   "RECURSIVE",
106:   "REMAP"
107: };

109: typedef struct {
110:   PetscSection section;   /* Sizes for each partition */
111:   IS           partition; /* Points in each partition */
112:   PetscBool    random;    /* Flag for a random partition */
113: } PetscPartitioner_Shell;

115: typedef struct {
116:   PetscInt dummy;
117: } PetscPartitioner_Simple;

119: typedef struct {
120:   PetscInt dummy;
121: } PetscPartitioner_Gather;

123: /* Utility struct to store the contents of a Fluent file in memory */
124: typedef struct {
125:   int          index;    /* Type of section */
126:   unsigned int zoneID;
127:   unsigned int first;
128:   unsigned int last;
129:   int          type;
130:   int          nd;       /* Either ND or element-type */
131:   void        *data;
132: } FluentSection;

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

146: /* Point Numbering in Plex:

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

150:    First Stratum:  Cells [height 0]
151:    Second Stratum: Vertices [depth 0]
152:    Third Stratum:  Faces [height 1]
153:    Fourth Stratum: Edges [depth 1]

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

164:   PetscSection         coneSection;       /* Layout of cones (inedges for DAG) */
165:   PetscInt             maxConeSize;       /* Cached for fast lookup */
166:   PetscInt            *cones;             /* Cone for each point */
167:   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 */
168:   PetscSection         supportSection;    /* Layout of cones (inedges for DAG) */
169:   PetscInt             maxSupportSize;    /* Cached for fast lookup */
170:   PetscInt            *supports;          /* Cone for each point */
171:   PetscBool            refinementUniform; /* Flag for uniform cell refinement */
172:   PetscReal            refinementLimit;   /* Maximum volume for refined cell */
173:   PetscErrorCode     (*refinementFunc)(const PetscReal [], PetscReal *); /* Function giving the maximum volume for refined cell */
174:   PetscInt             hybridPointMax[8]; /* Allow segregation of some points, each dimension has a divider (used in VTK output and refinement) */
175:   PetscInt             overlap;           /* Overlap of the partitions as passed to DMPlexDistribute() or DMPlexDistributeOverlap() */
176:   PetscInt             ghostCellStart;    /* The first ghost cell (for FV BC) or -1 */
177:   DMPlexInterpolatedFlag interpolated;
178:   DMPlexInterpolatedFlag interpolatedCollective;

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

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

185:   /* Generation */
186:   char                *tetgenOpts;
187:   char                *triangleOpts;
188:   PetscPartitioner     partitioner;
189:   PetscBool            partitionBalance;  /* Evenly divide partition overlap when distributing */
190:   PetscBool            remeshBd;

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

195:   /* Labels and numbering */
196:   PetscObjectState     depthState;        /* State of depth label, so that we can determine if a user changes it */
197:   IS                   globalVertexNumbers;
198:   IS                   globalCellNumbers;

200:   /* Constraints */
201:   PetscSection         anchorSection;      /* maps constrained points to anchor points */
202:   IS                   anchorIS;           /* anchors indexed by the above section */
203:   PetscErrorCode     (*createanchors)(DM); /* automatically compute anchors (probably from tree constraints) */
204:   PetscErrorCode     (*computeanchormatrix)(DM,PetscSection,PetscSection,Mat);

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

215:   /* MATIS support */
216:   PetscSection         subdomainSection;

218:   /* Adjacency */
219:   PetscBool            useAnchors;        /* Replace constrained points with their anchors in adjacency lists */
220:   PetscErrorCode      (*useradjacency)(DM,PetscInt,PetscInt*,PetscInt[],void*); /* User callback for adjacency */
221:   void                *useradjacencyctx;  /* User context for callback */

223:   /* Projection */
224:   PetscInt             maxProjectionHeight; /* maximum height of cells used in DMPlexProject functions */

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

230:   /* Geometry */
231:   PetscReal            minradius;         /* Minimum distance from cell centroid to face */
232:   PetscBool            useHashLocation;   /* Use grid hashing for point location */
233:   PetscGridHash        lbox;              /* Local box for searching */

235:   /* Debugging */
236:   PetscBool            printSetValues;
237:   PetscInt             printFEM;
238:   PetscInt             printL2;
239:   PetscReal            printTol;
240: } DM_Plex;

242: PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll_VTU(DM,PetscViewer);
243: PETSC_EXTERN PetscErrorCode VecView_Plex_Local(Vec,PetscViewer);
244: PETSC_EXTERN PetscErrorCode VecView_Plex_Native(Vec,PetscViewer);
245: PETSC_EXTERN PetscErrorCode VecView_Plex(Vec,PetscViewer);
246: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Local(Vec,PetscViewer);
247: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Native(Vec,PetscViewer);
248: PETSC_EXTERN PetscErrorCode VecLoad_Plex(Vec,PetscViewer);
249: PETSC_INTERN PetscErrorCode DMPlexGetFieldType_Internal(DM, PetscSection, PetscInt, PetscInt *, PetscInt *, PetscViewerVTKFieldType *);
250: PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM,PetscViewer);
251: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject,PetscViewer);
252: #if defined(PETSC_HAVE_HDF5)
253: PETSC_EXTERN PetscErrorCode VecView_Plex_Local_HDF5(Vec, PetscViewer);
254: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5(Vec, PetscViewer);
255: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5(Vec, PetscViewer);
256: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5_Native(Vec, PetscViewer);
257: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5_Native(Vec, PetscViewer);
258: PETSC_EXTERN PetscErrorCode DMPlexView_HDF5(DM, PetscViewer);
259: PETSC_EXTERN PetscErrorCode DMPlexLoad_HDF5(DM, PetscViewer);
260: #endif

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

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

325: PETSC_INTERN PetscErrorCode DMPlexCreateCellNumbering_Internal(DM, PetscBool, IS *);
326: PETSC_INTERN PetscErrorCode DMPlexCreateVertexNumbering_Internal(DM, PetscBool, IS *);
327: PETSC_INTERN PetscErrorCode DMPlexCreateNumbering_Internal(DM, PetscInt, PetscInt, PetscInt, PetscInt *, PetscSF, IS *);
328: PETSC_INTERN PetscErrorCode DMPlexRefine_Internal(DM, DMLabel, DM *);
329: PETSC_INTERN PetscErrorCode DMPlexCoarsen_Internal(DM, DMLabel, DM *);
330: PETSC_INTERN PetscErrorCode DMCreateMatrix_Plex(DM, Mat*);

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

334: /* invert dihedral symmetry: return a^-1,
335:  * using the representation described in
336:  * DMPlexGetConeOrientation() */
337: PETSC_STATIC_INLINE PetscInt DihedralInvert(PetscInt N, PetscInt a)
338: {
339:   return (a <= 0) ? a : (N - a);
340: }

342: /* invert dihedral symmetry: return b * a,
343:  * using the representation described in
344:  * DMPlexGetConeOrientation() */
345: PETSC_STATIC_INLINE PetscInt DihedralCompose(PetscInt N, PetscInt a, PetscInt b)
346: {
347:   if (!N) return 0;
348:   return  (a >= 0) ?
349:          ((b >= 0) ? ((a + b) % N) : -(((a - b - 1) % N) + 1)) :
350:          ((b >= 0) ? -(((N - b - a - 1) % N) + 1) : ((N + b - a) % N));
351: }

353: /* swap dihedral symmetries: return b * a^-1,
354:  * using the representation described in
355:  * DMPlexGetConeOrientation() */
356: PETSC_STATIC_INLINE PetscInt DihedralSwap(PetscInt N, PetscInt a, PetscInt b)
357: {
358:   return DihedralCompose(N,DihedralInvert(N,a),b);
359: }

361: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Internal(DM, IS , PetscReal, Vec, Vec, PetscReal, Vec, void *);
362: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Internal(DM, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
363: PETSC_EXTERN PetscErrorCode DMPlexReconstructGradients_Internal(DM, PetscFV, PetscInt, PetscInt, Vec, Vec, Vec, Vec);

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

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

480: PETSC_STATIC_INLINE void DMPlex_MatMultLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
481: {
482:   PetscInt j;
483:   for (j = 0; j < m; ++j) {
484:     DMPlex_MultTranspose2D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
485:   }
486:   (void)PetscLogFlops(8.0*m);
487: }
488: PETSC_STATIC_INLINE void DMPlex_MatMultLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
489: {
490:   PetscInt j;
491:   for (j = 0; j < m; ++j) {
492:     DMPlex_MultTranspose3D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
493:   }
494:   (void)PetscLogFlops(8.0*m);
495: }
496: PETSC_STATIC_INLINE void DMPlex_MatMultTransposeLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
497: {
498:   PetscInt j;
499:   for (j = 0; j < m; ++j) {
500:     DMPlex_Mult2D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
501:   }
502:   (void)PetscLogFlops(8.0*m);
503: }
504: PETSC_STATIC_INLINE void DMPlex_MatMultTransposeLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
505: {
506:   PetscInt j;
507:   for (j = 0; j < m; ++j) {
508:     DMPlex_Mult3D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
509:   }
510:   (void)PetscLogFlops(8.0*m);
511: }

513: PETSC_STATIC_INLINE void DMPlex_Transpose2D_Internal(PetscScalar A[])
514: {
515:   PetscScalar tmp;
516:   tmp = A[1]; A[1] = A[2]; A[2] = tmp;
517: }
518: PETSC_STATIC_INLINE void DMPlex_Transpose3D_Internal(PetscScalar A[])
519: {
520:   PetscScalar tmp;
521:   tmp = A[1]; A[1] = A[3]; A[3] = tmp;
522:   tmp = A[2]; A[2] = A[6]; A[6] = tmp;
523:   tmp = A[5]; A[5] = A[7]; A[7] = tmp;
524: }

526: PETSC_STATIC_INLINE void DMPlex_Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
527: {
528:   const PetscReal invDet = 1.0/detJ;

530:   invJ[0] =  invDet*J[3];
531:   invJ[1] = -invDet*J[1];
532:   invJ[2] = -invDet*J[2];
533:   invJ[3] =  invDet*J[0];
534:   (void)PetscLogFlops(5.0);
535: }

537: PETSC_STATIC_INLINE void DMPlex_Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
538: {
539:   const PetscReal invDet = 1.0/detJ;

541:   invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
542:   invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
543:   invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
544:   invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
545:   invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
546:   invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
547:   invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
548:   invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
549:   invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
550:   (void)PetscLogFlops(37.0);
551: }

553: PETSC_STATIC_INLINE void DMPlex_Det2D_Internal(PetscReal *detJ, const PetscReal J[])
554: {
555:   *detJ = J[0]*J[3] - J[1]*J[2];
556:   (void)PetscLogFlops(3.0);
557: }

559: PETSC_STATIC_INLINE void DMPlex_Det3D_Internal(PetscReal *detJ, const PetscReal J[])
560: {
561:   *detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
562:            J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
563:            J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
564:   (void)PetscLogFlops(12.0);
565: }

567: PETSC_STATIC_INLINE void DMPlex_Det2D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
568: {
569:   *detJ = PetscRealPart(J[0])*PetscRealPart(J[3]) - PetscRealPart(J[1])*PetscRealPart(J[2]);
570:   (void)PetscLogFlops(3.0);
571: }

573: PETSC_STATIC_INLINE void DMPlex_Det3D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
574: {
575:   *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])) +
576:            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])) +
577:            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])));
578:   (void)PetscLogFlops(12.0);
579: }

581: 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];}

583: 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;}

585: 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;}

587: 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);}

589: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Translate_Private(PetscInt ornt, PetscInt *start, PetscBool *reverse)
590: {
592:   *reverse = (ornt < 0) ? PETSC_TRUE : PETSC_FALSE;
593:   *start = *reverse ? -(ornt+1) : ornt;
594:   return(0);
595: }

597: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Combine_Private(PetscInt coneSize, PetscInt origStart, PetscBool origReverse, PetscInt rotateStart, PetscBool rotateReverse, PetscInt *newStart, PetscBool *newReverse)
598: {
600:   *newReverse = (origReverse == rotateReverse) ? PETSC_FALSE : PETSC_TRUE;
601:   *newStart = rotateReverse ? (coneSize + rotateStart - origStart) : (coneSize + origStart - rotateStart);
602:   *newStart %= coneSize;
603:   return(0);
604: }

606: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_TranslateBack_Private(PetscInt coneSize, PetscInt start, PetscBool reverse, PetscInt *ornt)
607: {
609:   if (coneSize < 3) {
610:     /* edges just get flipped if start == 1 regardless direction */
611:     *ornt = start ? -2 : 0;
612:   } else {
613:     *ornt = reverse ? -(start+1) : start;
614:   }
615:   return(0);
616: }

618: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Permute_Private(PetscInt n, const PetscInt arr[], PetscInt start, PetscBool reverse, PetscInt newarr[])
619: {
620:   PetscInt i;

623:   if (reverse) {for (i=0; i<n; i++) newarr[i] = arr[(n+start-i)%n];}
624:   else         {for (i=0; i<n; i++) newarr[i] = arr[(start+i)%n];}
625:   return(0);
626: }

628: PETSC_INTERN PetscErrorCode DMPlexGetPointDualSpaceFEM(DM,PetscInt,PetscInt,PetscDualSpace *);
629: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPoint_Internal(PetscSection,PetscInt,PetscInt,PetscInt *,PetscBool,const PetscInt[],const PetscInt[],PetscInt[]);
630: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPointFields_Internal(PetscSection,PetscInt,PetscInt,PetscInt[],PetscBool,const PetscInt***,PetscInt,const PetscInt[],PetscInt[]);
631: PETSC_INTERN PetscErrorCode DMPlexGetCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);
632: PETSC_INTERN PetscErrorCode DMPlexRestoreCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);

634: PETSC_EXTERN PetscErrorCode DMSNESGetFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
635: PETSC_EXTERN PetscErrorCode DMSNESRestoreFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
636: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Patch_Internal(DM, PetscSection, IS, PetscReal, Vec, Vec, Vec, void *);
637: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Patch_Internal(DM, PetscSection, PetscSection, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
638: PETSC_INTERN PetscErrorCode DMCreateSubDomainDM_Plex(DM,DMLabel,PetscInt,IS*,DM*);
639: PETSC_INTERN PetscErrorCode DMPlexBasisTransformPoint_Internal(DM, DM, Vec, PetscInt, PetscBool[], PetscBool, PetscScalar *);
640: PETSC_EXTERN PetscErrorCode DMPlexBasisTransformPointTensor_Internal(DM, DM, Vec, PetscInt, PetscBool, PetscInt, PetscScalar *);
641: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApplyReal_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscReal *, PetscReal *, void *);
642: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApply_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscScalar *, PetscScalar *, void *);
643: PETSC_INTERN PetscErrorCode DMCreateNeumannOverlap_Plex(DM, IS*, Mat*, PetscErrorCode (**)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*), void **);

645: #endif /* _PLEXIMPL_H */