Actual source code: dmpleximpl.h

petsc-master 2019-08-18
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: typedef struct {
147:   PetscInt             refct;

149:   /* Sieve */
150:   PetscSection         coneSection;       /* Layout of cones (inedges for DAG) */
151:   PetscInt             maxConeSize;       /* Cached for fast lookup */
152:   PetscInt            *cones;             /* Cone for each point */
153:   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 */
154:   PetscSection         supportSection;    /* Layout of cones (inedges for DAG) */
155:   PetscInt             maxSupportSize;    /* Cached for fast lookup */
156:   PetscInt            *supports;          /* Cone for each point */
157:   PetscBool            refinementUniform; /* Flag for uniform cell refinement */
158:   PetscReal            refinementLimit;   /* Maximum volume for refined cell */
159:   PetscErrorCode     (*refinementFunc)(const PetscReal [], PetscReal *); /* Function giving the maximum volume for refined cell */
160:   PetscInt             hybridPointMax[8]; /* Allow segregation of some points, each dimension has a divider (used in VTK output and refinement) */

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

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

167:   /* Generation */
168:   char                *tetgenOpts;
169:   char                *triangleOpts;
170:   PetscPartitioner     partitioner;
171:   PetscBool            partitionBalance;  /* Evenly divide partition overlap when distributing */
172:   PetscBool            remeshBd;

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

177:   /* Labels and numbering */
178:   PetscObjectState     depthState;        /* State of depth label, so that we can determine if a user changes it */
179:   IS                   globalVertexNumbers;
180:   IS                   globalCellNumbers;

182:   /* Constraints */
183:   PetscSection         anchorSection;      /* maps constrained points to anchor points */
184:   IS                   anchorIS;           /* anchors indexed by the above section */
185:   PetscErrorCode     (*createanchors)(DM); /* automatically compute anchors (probably from tree constraints) */
186:   PetscErrorCode     (*computeanchormatrix)(DM,PetscSection,PetscSection,Mat);

188:   /* Tree: automatically construct constraints for hierarchically non-conforming meshes */
189:   PetscSection         parentSection;     /* dof == 1 if point has parent */
190:   PetscInt            *parents;           /* point to parent */
191:   PetscInt            *childIDs;          /* point to child ID */
192:   PetscSection         childSection;      /* inverse of parent section */
193:   PetscInt            *children;          /* point to children */
194:   DM                   referenceTree;     /* reference tree to which child ID's refer */
195:   PetscErrorCode      (*getchildsymmetry)(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*);

197:   /* MATIS support */
198:   PetscSection         subdomainSection;

200:   /* Adjacency */
201:   PetscBool            useAnchors;        /* Replace constrained points with their anchors in adjacency lists */
202:   PetscErrorCode      (*useradjacency)(DM,PetscInt,PetscInt*,PetscInt[],void*); /* User callback for adjacency */
203:   void                *useradjacencyctx;  /* User context for callback */

205:   /* Projection */
206:   PetscInt             maxProjectionHeight; /* maximum height of cells used in DMPlexProject functions */

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

212:   /* Geometry */
213:   PetscReal            minradius;         /* Minimum distance from cell centroid to face */
214:   PetscBool            useHashLocation;   /* Use grid hashing for point location */
215:   PetscGridHash        lbox;              /* Local box for searching */

217:   /* Debugging */
218:   PetscBool            printSetValues;
219:   PetscInt             printFEM;
220:   PetscInt             printL2;
221:   PetscReal            printTol;
222: } DM_Plex;

224: PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll_VTU(DM,PetscViewer);
225: PETSC_EXTERN PetscErrorCode VecView_Plex_Local(Vec,PetscViewer);
226: PETSC_EXTERN PetscErrorCode VecView_Plex_Native(Vec,PetscViewer);
227: PETSC_EXTERN PetscErrorCode VecView_Plex(Vec,PetscViewer);
228: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Local(Vec,PetscViewer);
229: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Native(Vec,PetscViewer);
230: PETSC_EXTERN PetscErrorCode VecLoad_Plex(Vec,PetscViewer);
231: PETSC_INTERN PetscErrorCode DMPlexGetFieldType_Internal(DM, PetscSection, PetscInt, PetscInt *, PetscInt *, PetscViewerVTKFieldType *);
232: PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM,PetscViewer);
233: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject,PetscViewer);
234: #if defined(PETSC_HAVE_HDF5)
235: PETSC_EXTERN PetscErrorCode VecView_Plex_Local_HDF5(Vec, PetscViewer);
236: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5(Vec, PetscViewer);
237: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5(Vec, PetscViewer);
238: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5_Native(Vec, PetscViewer);
239: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5_Native(Vec, PetscViewer);
240: PETSC_EXTERN PetscErrorCode DMPlexView_HDF5(DM, PetscViewer);
241: PETSC_EXTERN PetscErrorCode DMPlexLoad_HDF5(DM, PetscViewer);
242: #endif

244: PETSC_INTERN PetscErrorCode DMPlexClosurePoints_Private(DM,PetscInt,const PetscInt[],IS*);
245: PETSC_INTERN PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *, DM);
246: PETSC_INTERN PetscErrorCode DMCoarsen_Plex(DM, MPI_Comm, DM *);
247: PETSC_INTERN PetscErrorCode DMCoarsenHierarchy_Plex(DM, PetscInt, DM []);
248: PETSC_INTERN PetscErrorCode DMRefine_Plex(DM, MPI_Comm, DM *);
249: PETSC_INTERN PetscErrorCode DMRefineHierarchy_Plex(DM, PetscInt, DM []);
250: PETSC_INTERN PetscErrorCode DMAdaptLabel_Plex(DM, DMLabel, DM *);
251: PETSC_INTERN PetscErrorCode DMAdaptMetric_Plex(DM, Vec, DMLabel, DM *);
252: PETSC_INTERN PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
253: PETSC_INTERN PetscErrorCode DMProjectFunctionLocal_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
254: 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);
255: 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);
256: 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);
257: PETSC_INTERN PetscErrorCode DMComputeL2Diff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
258: PETSC_INTERN PetscErrorCode DMComputeL2GradientDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[], const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,const PetscReal [],PetscReal *);
259: PETSC_INTERN PetscErrorCode DMComputeL2FieldDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
260: PETSC_INTERN PetscErrorCode DMLocatePoints_Plex(DM, Vec, DMPointLocationType, PetscSF);

262: PETSC_INTERN PetscErrorCode DMPlexBuildFromCellList_Internal(DM, PetscInt, PetscInt, PetscInt, PetscInt, const int[], PetscBool);
263: PETSC_INTERN PetscErrorCode DMPlexBuildFromCellList_Parallel_Internal(DM, PetscInt, PetscInt, PetscInt, PetscInt, const int[], PetscBool, PetscSF *);
264: PETSC_INTERN PetscErrorCode DMPlexBuildCoordinates_Internal(DM, PetscInt, PetscInt, PetscInt, const double[]);
265: PETSC_INTERN PetscErrorCode DMPlexBuildCoordinates_Parallel_Internal(DM, PetscInt, PetscInt, PetscInt, PetscSF, const PetscReal[]);
266: PETSC_INTERN PetscErrorCode DMPlexLoadLabels_HDF5_Internal(DM, PetscViewer);
267: PETSC_INTERN PetscErrorCode DMPlexView_HDF5_Internal(DM, PetscViewer);
268: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Internal(DM, PetscViewer);
269: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Xdmf_Internal(DM, PetscViewer);
270: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Internal(Vec, PetscViewer);
271: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Native_Internal(Vec, PetscViewer);
272: PETSC_INTERN PetscErrorCode VecView_Plex_Local_HDF5_Internal(Vec, PetscViewer);
273: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Internal(Vec, PetscViewer);
274: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Native_Internal(Vec, PetscViewer);
275: /* TODO Make these INTERN */
276: PETSC_EXTERN PetscErrorCode DMPlexView_ExodusII_Internal(DM, int, PetscInt);
277: PETSC_EXTERN PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec, int, int);
278: PETSC_EXTERN PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec, int, int);
279: PETSC_EXTERN PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec, int, int);
280: PETSC_EXTERN PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec, int, int);
281: PETSC_INTERN PetscErrorCode DMPlexVTKGetCellType_Internal(DM,PetscInt,PetscInt,PetscInt*);
282: PETSC_INTERN PetscErrorCode DMPlexGetAdjacency_Internal(DM,PetscInt,PetscBool,PetscBool,PetscBool,PetscInt*,PetscInt*[]);
283: PETSC_INTERN PetscErrorCode DMPlexGetFaces_Internal(DM,PetscInt,PetscInt,PetscInt*,PetscInt*,const PetscInt*[]);
284: PETSC_INTERN PetscErrorCode DMPlexGetRawFaces_Internal(DM,PetscInt,PetscInt,const PetscInt[], PetscInt*,PetscInt*,const PetscInt*[]);
285: PETSC_INTERN PetscErrorCode DMPlexRestoreFaces_Internal(DM,PetscInt,PetscInt,PetscInt*,PetscInt*,const PetscInt*[]);
286: PETSC_INTERN PetscErrorCode DMPlexRefineUniform_Internal(DM,CellRefiner,DM*);
287: PETSC_INTERN PetscErrorCode DMPlexGetCellRefiner_Internal(DM,CellRefiner*);
288: PETSC_INTERN PetscErrorCode CellRefinerGetAffineTransforms_Internal(CellRefiner, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[]);
289: PETSC_INTERN PetscErrorCode CellRefinerGetAffineFaceTransforms_Internal(CellRefiner, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[], PetscReal *[]);
290: PETSC_INTERN PetscErrorCode CellRefinerRestoreAffineTransforms_Internal(CellRefiner, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[]);
291: PETSC_INTERN PetscErrorCode CellRefinerInCellTest_Internal(CellRefiner, const PetscReal[], PetscBool *);
292: PETSC_INTERN PetscErrorCode DMPlexInvertCell_Internal(PetscInt, PetscInt, PetscInt[]);
293: PETSC_INTERN PetscErrorCode DMPlexVecSetFieldClosure_Internal(DM, PetscSection, Vec, PetscBool[], PetscInt, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
294: PETSC_INTERN PetscErrorCode DMPlexProjectConstraints_Internal(DM, Vec, Vec);
295: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_SetTree(DM, PetscSection, PetscInt[], PetscInt[]);
296: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_Union(DM,DM,const char *,DM*);
297: PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorTree(DM,DM,PetscSF,PetscInt *,Mat);
298: PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorTree(DM,DM,PetscSF,PetscInt *,Mat);
299: PETSC_EXTERN PetscErrorCode DMPlexAnchorsModifyMat(DM,PetscSection,PetscInt,PetscInt,const PetscInt[],const PetscInt ***,const PetscScalar[],PetscInt*,PetscInt*,PetscInt*[],PetscScalar*[],PetscInt[],PetscBool);
300: PETSC_EXTERN PetscErrorCode indicesPoint_private(PetscSection,PetscInt,PetscInt,PetscInt *,PetscBool,PetscInt,PetscInt []);
301: PETSC_EXTERN PetscErrorCode indicesPointFields_private(PetscSection,PetscInt,PetscInt,PetscInt [],PetscBool,PetscInt,PetscInt []);
302: PETSC_INTERN PetscErrorCode DMPlexLocatePoint_Internal(DM,PetscInt,const PetscScalar [],PetscInt,PetscInt *);
303: PETSC_EXTERN PetscErrorCode DMPlexOrientCell_Internal(DM,PetscInt,PetscInt,PetscBool);
304: PETSC_EXTERN PetscErrorCode DMPlexOrientInterface(DM);

306: PETSC_INTERN PetscErrorCode DMPlexCreateCellNumbering_Internal(DM, PetscBool, IS *);
307: PETSC_INTERN PetscErrorCode DMPlexCreateVertexNumbering_Internal(DM, PetscBool, IS *);
308: PETSC_INTERN PetscErrorCode DMPlexCreateNumbering_Internal(DM, PetscInt, PetscInt, PetscInt, PetscInt *, PetscSF, IS *);
309: PETSC_INTERN PetscErrorCode DMPlexRefine_Internal(DM, DMLabel, DM *);
310: PETSC_INTERN PetscErrorCode DMPlexCoarsen_Internal(DM, DMLabel, DM *);
311: PETSC_INTERN PetscErrorCode DMCreateMatrix_Plex(DM, Mat*);

313: /* invert dihedral symmetry: return a^-1,
314:  * using the representation described in
315:  * DMPlexGetConeOrientation() */
316: PETSC_STATIC_INLINE PetscInt DihedralInvert(PetscInt N, PetscInt a)
317: {
318:   return (a <= 0) ? a : (N - a);
319: }

321: /* invert dihedral symmetry: return b * a,
322:  * using the representation described in
323:  * DMPlexGetConeOrientation() */
324: PETSC_STATIC_INLINE PetscInt DihedralCompose(PetscInt N, PetscInt a, PetscInt b)
325: {
326:   if (!N) return 0;
327:   return  (a >= 0) ?
328:          ((b >= 0) ? ((a + b) % N) : -(((a - b - 1) % N) + 1)) :
329:          ((b >= 0) ? -(((N - b - a - 1) % N) + 1) : ((N + b - a) % N));
330: }

332: /* swap dihedral symmetries: return b * a^-1,
333:  * using the representation described in
334:  * DMPlexGetConeOrientation() */
335: PETSC_STATIC_INLINE PetscInt DihedralSwap(PetscInt N, PetscInt a, PetscInt b)
336: {
337:   return DihedralCompose(N,DihedralInvert(N,a),b);
338: }

340: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Internal(DM, IS , PetscReal, Vec, Vec, PetscReal, Vec, void *);
341: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Internal(DM, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
342: PETSC_EXTERN PetscErrorCode DMPlexReconstructGradients_Internal(DM, PetscFV, PetscInt, PetscInt, Vec, Vec, Vec, Vec);

344: /* Matvec with A in row-major storage, x and y can be aliased */
345: PETSC_STATIC_INLINE void DMPlex_Mult2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
346: {
347:   PetscScalar z[2];
348:   z[0] = x[0]; z[1] = x[ldx];
349:   y[0]   = A[0]*z[0] + A[1]*z[1];
350:   y[ldx] = A[2]*z[0] + A[3]*z[1];
351:   (void)PetscLogFlops(6.0);
352: }
353: PETSC_STATIC_INLINE void DMPlex_Mult3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
354: {
355:   PetscScalar z[3];
356:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
357:   y[0]     = A[0]*z[0] + A[1]*z[1] + A[2]*z[2];
358:   y[ldx]   = A[3]*z[0] + A[4]*z[1] + A[5]*z[2];
359:   y[ldx*2] = A[6]*z[0] + A[7]*z[1] + A[8]*z[2];
360:   (void)PetscLogFlops(15.0);
361: }
362: PETSC_STATIC_INLINE void DMPlex_MultTranspose2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
363: {
364:   PetscScalar z[2];
365:   z[0] = x[0]; z[1] = x[ldx];
366:   y[0]   = A[0]*z[0] + A[2]*z[1];
367:   y[ldx] = A[1]*z[0] + A[3]*z[1];
368:   (void)PetscLogFlops(6.0);
369: }
370: PETSC_STATIC_INLINE void DMPlex_MultTranspose3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
371: {
372:   PetscScalar z[3];
373:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
374:   y[0]     = A[0]*z[0] + A[3]*z[1] + A[6]*z[2];
375:   y[ldx]   = A[1]*z[0] + A[4]*z[1] + A[7]*z[2];
376:   y[ldx*2] = A[2]*z[0] + A[5]*z[1] + A[8]*z[2];
377:   (void)PetscLogFlops(15.0);
378: }
379: PETSC_STATIC_INLINE void DMPlex_Mult2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
380: {
381:   PetscScalar z[2];
382:   z[0] = x[0]; z[1] = x[ldx];
383:   y[0]   = A[0]*z[0] + A[1]*z[1];
384:   y[ldx] = A[2]*z[0] + A[3]*z[1];
385:   (void)PetscLogFlops(6.0);
386: }
387: PETSC_STATIC_INLINE void DMPlex_Mult3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
388: {
389:   PetscScalar z[3];
390:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
391:   y[0]     = A[0]*z[0] + A[1]*z[1] + A[2]*z[2];
392:   y[ldx]   = A[3]*z[0] + A[4]*z[1] + A[5]*z[2];
393:   y[ldx*2] = A[6]*z[0] + A[7]*z[1] + A[8]*z[2];
394:   (void)PetscLogFlops(15.0);
395: }
396: PETSC_STATIC_INLINE void DMPlex_MultTranspose2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
397: {
398:   PetscScalar z[2];
399:   z[0] = x[0]; z[1] = x[ldx];
400:   y[0]   = A[0]*z[0] + A[2]*z[1];
401:   y[ldx] = A[1]*z[0] + A[3]*z[1];
402:   (void)PetscLogFlops(6.0);
403: }
404: PETSC_STATIC_INLINE void DMPlex_MultTranspose3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
405: {
406:   PetscScalar z[3];
407:   z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
408:   y[0]     = A[0]*z[0] + A[3]*z[1] + A[6]*z[2];
409:   y[ldx]   = A[1]*z[0] + A[4]*z[1] + A[7]*z[2];
410:   y[ldx*2] = A[2]*z[0] + A[5]*z[1] + A[8]*z[2];
411:   (void)PetscLogFlops(15.0);
412: }

414: PETSC_STATIC_INLINE void DMPlex_MatMult2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
415: {
416:   PetscInt j;
417:   for (j = 0; j < n; ++j) {
418:     PetscScalar z[2];
419:     z[0] = B[0+j]; z[1] = B[1*ldb+j];
420:     DMPlex_Mult2D_Internal(A, 1, z, z);
421:     C[0+j] = z[0]; C[1*ldb+j] = z[1];
422:   }
423:   (void)PetscLogFlops(8.0*n);
424: }
425: PETSC_STATIC_INLINE void DMPlex_MatMult3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
426: {
427:   PetscInt j;
428:   for (j = 0; j < n; ++j) {
429:     PetscScalar z[3];
430:     z[0] = B[0+j]; z[1] = B[1*ldb+j]; z[2] = B[2*ldb+j];
431:     DMPlex_Mult3D_Internal(A, 1, z, z);
432:     C[0+j] = z[0]; C[1*ldb+j] = z[1]; C[2*ldb+j] = z[2];
433:   }
434:   (void)PetscLogFlops(8.0*n);
435: }
436: PETSC_STATIC_INLINE void DMPlex_MatMultTranspose2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
437: {
438:   PetscInt j;
439:   for (j = 0; j < n; ++j) {
440:     PetscScalar z[2];
441:     z[0] = B[0+j]; z[1] = B[1*ldb+j];
442:     DMPlex_MultTranspose2D_Internal(A, 1, z, z);
443:     C[0+j] = z[0]; C[1*ldb+j] = z[1];
444:   }
445:   (void)PetscLogFlops(8.0*n);
446: }
447: PETSC_STATIC_INLINE void DMPlex_MatMultTranspose3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
448: {
449:   PetscInt j;
450:   for (j = 0; j < n; ++j) {
451:     PetscScalar z[3];
452:     z[0] = B[0+j]; z[1] = B[1*ldb+j]; z[2] = B[2*ldb+j];
453:     DMPlex_MultTranspose3D_Internal(A, 1, z, z);
454:     C[0+j] = z[0]; C[1*ldb+j] = z[1]; C[2*ldb+j] = z[2];
455:   }
456:   (void)PetscLogFlops(8.0*n);
457: }

459: PETSC_STATIC_INLINE void DMPlex_MatMultLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
460: {
461:   PetscInt j;
462:   for (j = 0; j < m; ++j) {
463:     DMPlex_MultTranspose2D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
464:   }
465:   (void)PetscLogFlops(8.0*m);
466: }
467: PETSC_STATIC_INLINE void DMPlex_MatMultLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
468: {
469:   PetscInt j;
470:   for (j = 0; j < m; ++j) {
471:     DMPlex_MultTranspose3D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
472:   }
473:   (void)PetscLogFlops(8.0*m);
474: }
475: PETSC_STATIC_INLINE void DMPlex_MatMultTransposeLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
476: {
477:   PetscInt j;
478:   for (j = 0; j < m; ++j) {
479:     DMPlex_Mult2D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
480:   }
481:   (void)PetscLogFlops(8.0*m);
482: }
483: PETSC_STATIC_INLINE void DMPlex_MatMultTransposeLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
484: {
485:   PetscInt j;
486:   for (j = 0; j < m; ++j) {
487:     DMPlex_Mult3D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
488:   }
489:   (void)PetscLogFlops(8.0*m);
490: }

492: PETSC_STATIC_INLINE void DMPlex_Transpose2D_Internal(PetscScalar A[])
493: {
494:   PetscScalar tmp;
495:   tmp = A[1]; A[1] = A[2]; A[2] = tmp;
496: }
497: PETSC_STATIC_INLINE void DMPlex_Transpose3D_Internal(PetscScalar A[])
498: {
499:   PetscScalar tmp;
500:   tmp = A[1]; A[1] = A[3]; A[3] = tmp;
501:   tmp = A[2]; A[2] = A[6]; A[6] = tmp;
502:   tmp = A[5]; A[5] = A[7]; A[7] = tmp;
503: }

505: PETSC_STATIC_INLINE void DMPlex_Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
506: {
507:   const PetscReal invDet = 1.0/detJ;

509:   invJ[0] =  invDet*J[3];
510:   invJ[1] = -invDet*J[1];
511:   invJ[2] = -invDet*J[2];
512:   invJ[3] =  invDet*J[0];
513:   (void)PetscLogFlops(5.0);
514: }

516: PETSC_STATIC_INLINE void DMPlex_Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
517: {
518:   const PetscReal invDet = 1.0/detJ;

520:   invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
521:   invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
522:   invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
523:   invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
524:   invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
525:   invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
526:   invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
527:   invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
528:   invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
529:   (void)PetscLogFlops(37.0);
530: }

532: PETSC_STATIC_INLINE void DMPlex_Det2D_Internal(PetscReal *detJ, const PetscReal J[])
533: {
534:   *detJ = J[0]*J[3] - J[1]*J[2];
535:   (void)PetscLogFlops(3.0);
536: }

538: PETSC_STATIC_INLINE void DMPlex_Det3D_Internal(PetscReal *detJ, const PetscReal J[])
539: {
540:   *detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
541:            J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
542:            J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
543:   (void)PetscLogFlops(12.0);
544: }

546: PETSC_STATIC_INLINE void DMPlex_Det2D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
547: {
548:   *detJ = PetscRealPart(J[0])*PetscRealPart(J[3]) - PetscRealPart(J[1])*PetscRealPart(J[2]);
549:   (void)PetscLogFlops(3.0);
550: }

552: PETSC_STATIC_INLINE void DMPlex_Det3D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
553: {
554:   *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])) +
555:            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])) +
556:            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])));
557:   (void)PetscLogFlops(12.0);
558: }

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

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

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

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

568: /* Compare cones of the master and slave face (with the same cone points modulo order), and return relative orientation of the slave. */
569: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Orient_Private(PetscInt coneSize, PetscInt masterConeSize, const PetscInt masterCone[], const PetscInt slaveCone[], PetscInt *start, PetscBool *reverse)
570: {
571:   PetscInt        i;

574:   *start = 0;
575:   for (i=0; i<coneSize; i++) {
576:     if (slaveCone[i] == masterCone[0]) {
577:       *start = i;
578:       break;
579:     }
580:   }
581:   if (PetscUnlikely(i==coneSize)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "starting point of master cone not found in slave cone");
582:   *reverse = PETSC_FALSE;
583:   for (i=0; i<masterConeSize; i++) {if (slaveCone[((*start)+i)%coneSize] != masterCone[i]) break;}
584:   if (i == masterConeSize) return(0);
585:   *reverse = PETSC_TRUE;
586:   for (i=0; i<masterConeSize; i++) {if (slaveCone[(coneSize+(*start)-i)%coneSize] != masterCone[i]) break;}
587:   if (i < masterConeSize) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "master and slave cone have non-conforming order of points");
588:   return(0);
589: }

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

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

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

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

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

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

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

646: #endif /* _PLEXIMPL_H */