Actual source code: dmpleximpl.h
petsc-3.13.0 2020-03-29
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: typedef struct _PetscPartitionerOps *PetscPartitionerOps;
48: struct _PetscPartitionerOps {
49: PetscErrorCode (*setfromoptions)(PetscOptionItems*,PetscPartitioner);
50: PetscErrorCode (*setup)(PetscPartitioner);
51: PetscErrorCode (*view)(PetscPartitioner,PetscViewer);
52: PetscErrorCode (*destroy)(PetscPartitioner);
53: PetscErrorCode (*partition)(PetscPartitioner, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscSection, PetscSection, PetscSection, IS *);
54: };
56: struct _p_PetscPartitioner {
57: PETSCHEADER(struct _PetscPartitionerOps);
58: void *data; /* Implementation object */
59: PetscInt height; /* Height of points to partition into non-overlapping subsets */
60: PetscInt edgeCut; /* The number of edge cut by the partition */
61: PetscReal balance; /* The maximum partition size divided by the minimum size */
62: PetscViewer viewer;
63: PetscViewer viewerGraph;
64: PetscBool viewGraph;
65: PetscBool noGraph; /* if true, the partitioner does not need the connectivity graph, only the number of local vertices */
66: PetscBool usevwgt; /* if true, the partitioner looks at the local section vertSection to weight the vertices of the graph */
67: };
69: typedef struct {
70: PetscInt dummy;
71: } PetscPartitioner_Chaco;
73: typedef struct {
74: MPI_Comm pcomm;
75: PetscInt ptype;
76: PetscReal imbalanceRatio;
77: PetscInt debugFlag;
78: PetscInt randomSeed;
79: } PetscPartitioner_ParMetis;
81: typedef struct {
82: MPI_Comm pcomm;
83: PetscInt strategy;
84: PetscReal imbalance;
85: } PetscPartitioner_PTScotch;
87: static const char *const
88: PTScotchStrategyList[] = {
89: "DEFAULT",
90: "QUALITY",
91: "SPEED",
92: "BALANCE",
93: "SAFETY",
94: "SCALABILITY",
95: "RECURSIVE",
96: "REMAP"
97: };
99: typedef struct {
100: PetscSection section; /* Sizes for each partition */
101: IS partition; /* Points in each partition */
102: PetscBool random; /* Flag for a random partition */
103: } PetscPartitioner_Shell;
105: typedef struct {
106: PetscInt dummy;
107: } PetscPartitioner_Simple;
109: typedef struct {
110: PetscInt dummy;
111: } PetscPartitioner_Gather;
113: typedef struct _DMPlexCellRefinerOps *DMPlexCellRefinerOps;
114: struct _DMPlexCellRefinerOps {
115: PetscErrorCode (*refine)(DMPlexCellRefiner, DMPolytopeType, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]);
116: PetscErrorCode (*mapsubcells)(DMPlexCellRefiner, DMPolytopeType, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *);
117: PetscErrorCode (*getaffinetransforms)(DMPlexCellRefiner, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[]);
118: PetscErrorCode (*getaffinefacetransforms)(DMPlexCellRefiner, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[], PetscReal *[]);
119: PetscErrorCode (*getcellvertices)(DMPlexCellRefiner, DMPolytopeType, PetscInt *, PetscReal *[]);
120: PetscErrorCode (*getsubcellvertices)(DMPlexCellRefiner, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt *, PetscInt *[]);
121: };
123: struct _p_DMPlexCellRefiner {
124: PETSCHEADER(struct _DMPlexCellRefinerOps);
125: DM dm; /* The original DM */
126: PetscBool setupcalled;
127: DMPlexCellRefinerType type;
128: PetscInt *ctOrder; /* [i] = ct: An array with cell types in depth order */
129: PetscInt *ctOrderInv; /* [ct] = i: An array with the ordinal numbers for each cell type */
130: PetscInt *ctStart; /* The number for the first cell of each polytope type in the original mesh, indexed by cell type */
131: PetscInt *ctStartNew; /* The number for the first cell of each polytope type in the new mesh, indexed by cell type */
132: PetscInt *offset; /* [ct][ctNew]: The offset in the new point numbering of a point of type ctNew produced from an old point of type ct */
133: PetscFE *coordFE; /* Finite element for each cell type, used for localized coordinate interpolation */
134: PetscFEGeom **refGeom; /* Geometry of the reference cell for each cell type */
135: };
137: /* Utility struct to store the contents of a Fluent file in memory */
138: typedef struct {
139: int index; /* Type of section */
140: unsigned int zoneID;
141: unsigned int first;
142: unsigned int last;
143: int type;
144: int nd; /* Either ND or element-type */
145: void *data;
146: } FluentSection;
148: struct _PetscGridHash {
149: PetscInt dim;
150: PetscReal lower[3]; /* The lower-left corner */
151: PetscReal upper[3]; /* The upper-right corner */
152: PetscReal extent[3]; /* The box size */
153: PetscReal h[3]; /* The subbox size */
154: PetscInt n[3]; /* The number of subboxes */
155: PetscSection cellSection; /* Offsets for cells in each subbox*/
156: IS cells; /* List of cells in each subbox */
157: DMLabel cellsSparse; /* Sparse storage for cell map */
158: };
160: /* Point Numbering in Plex:
162: Points are numbered contiguously by stratum. Strate are organized as follows:
164: First Stratum: Cells [height 0]
165: Second Stratum: Vertices [depth 0]
166: Third Stratum: Faces [height 1]
167: Fourth Stratum: Edges [depth 1]
169: We do this so that the numbering of a cell-vertex mesh does not change after interpolation. Within a given stratum,
170: we allow additional segregation of by cell type.
171: */
172: typedef struct {
173: PetscInt refct;
175: PetscSection coneSection; /* Layout of cones (inedges for DAG) */
176: PetscInt maxConeSize; /* Cached for fast lookup */
177: PetscInt *cones; /* Cone for each point */
178: 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 */
179: PetscSection supportSection; /* Layout of cones (inedges for DAG) */
180: PetscInt maxSupportSize; /* Cached for fast lookup */
181: PetscInt *supports; /* Cone for each point */
182: PetscBool refinementUniform; /* Flag for uniform cell refinement */
183: PetscReal refinementLimit; /* Maximum volume for refined cell */
184: PetscErrorCode (*refinementFunc)(const PetscReal [], PetscReal *); /* Function giving the maximum volume for refined cell */
185: PetscInt overlap; /* Overlap of the partitions as passed to DMPlexDistribute() or DMPlexDistributeOverlap() */
186: DMPlexInterpolatedFlag interpolated;
187: DMPlexInterpolatedFlag interpolatedCollective;
189: PetscInt *facesTmp; /* Work space for faces operation */
191: /* Hierarchy */
192: DMPlexCellRefinerType cellRefiner; /* Strategy for refining cells */
193: PetscBool regularRefinement; /* This flag signals that we are a regular refinement of coarseMesh */
195: /* Generation */
196: char *tetgenOpts;
197: char *triangleOpts;
198: PetscPartitioner partitioner;
199: PetscBool partitionBalance; /* Evenly divide partition overlap when distributing */
200: PetscBool remeshBd;
202: /* Submesh */
203: DMLabel subpointMap; /* Label each original mesh point in the submesh with its depth, subpoint are the implicit numbering */
205: /* Labels and numbering */
206: PetscObjectState depthState; /* State of depth label, so that we can determine if a user changes it */
207: PetscObjectState celltypeState; /* State of celltype label, so that we can determine if a user changes it */
208: IS globalVertexNumbers;
209: IS globalCellNumbers;
211: /* Constraints */
212: PetscSection anchorSection; /* maps constrained points to anchor points */
213: IS anchorIS; /* anchors indexed by the above section */
214: PetscErrorCode (*createanchors)(DM); /* automatically compute anchors (probably from tree constraints) */
215: PetscErrorCode (*computeanchormatrix)(DM,PetscSection,PetscSection,Mat);
217: /* Tree: automatically construct constraints for hierarchically non-conforming meshes */
218: PetscSection parentSection; /* dof == 1 if point has parent */
219: PetscInt *parents; /* point to parent */
220: PetscInt *childIDs; /* point to child ID */
221: PetscSection childSection; /* inverse of parent section */
222: PetscInt *children; /* point to children */
223: DM referenceTree; /* reference tree to which child ID's refer */
224: PetscErrorCode (*getchildsymmetry)(DM,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*);
226: /* MATIS support */
227: PetscSection subdomainSection;
229: /* Adjacency */
230: PetscBool useAnchors; /* Replace constrained points with their anchors in adjacency lists */
231: PetscErrorCode (*useradjacency)(DM,PetscInt,PetscInt*,PetscInt[],void*); /* User callback for adjacency */
232: void *useradjacencyctx; /* User context for callback */
234: /* Projection */
235: PetscInt maxProjectionHeight; /* maximum height of cells used in DMPlexProject functions */
237: /* Output */
238: PetscInt vtkCellHeight; /* The height of cells for output, default is 0 */
239: PetscReal scale[NUM_PETSC_UNITS]; /* The scale for each SI unit */
241: /* Geometry */
242: PetscReal minradius; /* Minimum distance from cell centroid to face */
243: PetscBool useHashLocation; /* Use grid hashing for point location */
244: PetscGridHash lbox; /* Local box for searching */
246: /* Debugging */
247: PetscBool printSetValues;
248: PetscInt printFEM;
249: PetscInt printL2;
250: PetscReal printTol;
251: } DM_Plex;
253: PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll_VTU(DM,PetscViewer);
254: PETSC_EXTERN PetscErrorCode VecView_Plex_Local(Vec,PetscViewer);
255: PETSC_EXTERN PetscErrorCode VecView_Plex_Native(Vec,PetscViewer);
256: PETSC_EXTERN PetscErrorCode VecView_Plex(Vec,PetscViewer);
257: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Local(Vec,PetscViewer);
258: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Native(Vec,PetscViewer);
259: PETSC_EXTERN PetscErrorCode VecLoad_Plex(Vec,PetscViewer);
260: PETSC_INTERN PetscErrorCode DMPlexGetFieldType_Internal(DM, PetscSection, PetscInt, PetscInt *, PetscInt *, PetscViewerVTKFieldType *);
261: PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM,PetscViewer);
262: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject,PetscViewer);
263: #if defined(PETSC_HAVE_HDF5)
264: PETSC_EXTERN PetscErrorCode VecView_Plex_Local_HDF5(Vec, PetscViewer);
265: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5(Vec, PetscViewer);
266: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5(Vec, PetscViewer);
267: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5_Native(Vec, PetscViewer);
268: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5_Native(Vec, PetscViewer);
269: PETSC_EXTERN PetscErrorCode DMPlexView_HDF5(DM, PetscViewer);
270: PETSC_EXTERN PetscErrorCode DMPlexLoad_HDF5(DM, PetscViewer);
271: #endif
273: PETSC_INTERN PetscErrorCode DMPlexClosurePoints_Private(DM,PetscInt,const PetscInt[],IS*);
274: PETSC_INTERN PetscErrorCode DMSetFromOptions_NonRefinement_Plex(PetscOptionItems *, DM);
275: PETSC_INTERN PetscErrorCode DMCoarsen_Plex(DM, MPI_Comm, DM *);
276: PETSC_INTERN PetscErrorCode DMCoarsenHierarchy_Plex(DM, PetscInt, DM []);
277: PETSC_INTERN PetscErrorCode DMRefine_Plex(DM, MPI_Comm, DM *);
278: PETSC_INTERN PetscErrorCode DMRefineHierarchy_Plex(DM, PetscInt, DM []);
279: PETSC_INTERN PetscErrorCode DMAdaptLabel_Plex(DM, DMLabel, DM *);
280: PETSC_INTERN PetscErrorCode DMAdaptMetric_Plex(DM, Vec, DMLabel, DM *);
281: PETSC_INTERN PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
282: PETSC_INTERN PetscErrorCode DMProjectFunctionLocal_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
283: 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);
284: 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);
285: 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);
286: PETSC_INTERN PetscErrorCode DMComputeL2Diff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
287: PETSC_INTERN PetscErrorCode DMComputeL2GradientDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[], const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,const PetscReal [],PetscReal *);
288: PETSC_INTERN PetscErrorCode DMComputeL2FieldDiff_Plex(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
289: PETSC_INTERN PetscErrorCode DMLocatePoints_Plex(DM, Vec, DMPointLocationType, PetscSF);
291: PETSC_INTERN PetscErrorCode DMPlexBuildFromCellList_Internal(DM, PetscInt, PetscInt, PetscInt, PetscInt, const int[], PetscBool);
292: PETSC_INTERN PetscErrorCode DMPlexBuildFromCellList_Parallel_Internal(DM, PetscInt, PetscInt, PetscInt, PetscInt, const int[], PetscBool, PetscSF *);
293: PETSC_INTERN PetscErrorCode DMPlexBuildCoordinates_Internal(DM, PetscInt, PetscInt, PetscInt, const double[]);
294: PETSC_INTERN PetscErrorCode DMPlexBuildCoordinates_Parallel_Internal(DM, PetscInt, PetscInt, PetscInt, PetscSF, const PetscReal[]);
295: PETSC_INTERN PetscErrorCode DMPlexLoadLabels_HDF5_Internal(DM, PetscViewer);
296: PETSC_INTERN PetscErrorCode DMPlexView_HDF5_Internal(DM, PetscViewer);
297: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Internal(DM, PetscViewer);
298: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Xdmf_Internal(DM, PetscViewer);
299: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Internal(Vec, PetscViewer);
300: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Native_Internal(Vec, PetscViewer);
301: PETSC_INTERN PetscErrorCode VecView_Plex_Local_HDF5_Internal(Vec, PetscViewer);
302: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Internal(Vec, PetscViewer);
303: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Native_Internal(Vec, PetscViewer);
304: /* TODO Make these INTERN */
305: PETSC_EXTERN PetscErrorCode DMPlexView_ExodusII_Internal(DM, int, PetscInt);
306: PETSC_EXTERN PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec, int, int);
307: PETSC_EXTERN PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec, int, int);
308: PETSC_EXTERN PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec, int, int);
309: PETSC_EXTERN PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec, int, int);
310: PETSC_INTERN PetscErrorCode DMPlexVTKGetCellType_Internal(DM,PetscInt,PetscInt,PetscInt*);
311: PETSC_INTERN PetscErrorCode DMPlexGetAdjacency_Internal(DM,PetscInt,PetscBool,PetscBool,PetscBool,PetscInt*,PetscInt*[]);
312: PETSC_INTERN PetscErrorCode DMPlexGetRawFaces_Internal(DM,DMPolytopeType,const PetscInt[],PetscInt*,const DMPolytopeType*[],const PetscInt*[],const PetscInt*[]);
313: PETSC_INTERN PetscErrorCode DMPlexRestoreRawFaces_Internal(DM,DMPolytopeType,const PetscInt[],PetscInt*,const DMPolytopeType*[],const PetscInt*[],const PetscInt*[]);
314: PETSC_INTERN PetscErrorCode CellRefinerInCellTest_Internal(DMPolytopeType, const PetscReal[], PetscBool *);
315: PETSC_INTERN PetscErrorCode DMPlexComputeCellType_Internal(DM, PetscInt, PetscInt, DMPolytopeType *);
316: PETSC_INTERN PetscErrorCode DMPlexCreateCellTypeOrder_Internal(DMPolytopeType, PetscInt *[], PetscInt *[]);
317: PETSC_INTERN PetscErrorCode DMPlexVecSetFieldClosure_Internal(DM, PetscSection, Vec, PetscBool[], PetscInt, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
318: PETSC_INTERN PetscErrorCode DMPlexProjectConstraints_Internal(DM, Vec, Vec);
319: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_SetTree(DM, PetscSection, PetscInt[], PetscInt[]);
320: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_Union(DM,DM,const char *,DM*);
321: PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorTree(DM,DM,PetscSF,PetscInt *,Mat);
322: PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorTree(DM,DM,PetscSF,PetscInt *,Mat);
323: PETSC_EXTERN PetscErrorCode DMPlexAnchorsModifyMat(DM,PetscSection,PetscInt,PetscInt,const PetscInt[],const PetscInt ***,const PetscScalar[],PetscInt*,PetscInt*,PetscInt*[],PetscScalar*[],PetscInt[],PetscBool);
324: PETSC_EXTERN PetscErrorCode indicesPoint_private(PetscSection,PetscInt,PetscInt,PetscInt *,PetscBool,PetscInt,PetscInt []);
325: PETSC_EXTERN PetscErrorCode indicesPointFields_private(PetscSection,PetscInt,PetscInt,PetscInt [],PetscBool,PetscInt,PetscInt []);
326: PETSC_INTERN PetscErrorCode DMPlexLocatePoint_Internal(DM,PetscInt,const PetscScalar [],PetscInt,PetscInt *);
327: /* these two are PETSC_EXTERN just because of src/dm/impls/plex/tests/ex18.c */
328: PETSC_EXTERN PetscErrorCode DMPlexOrientCell_Internal(DM,PetscInt,PetscInt,PetscBool);
329: PETSC_EXTERN PetscErrorCode DMPlexOrientInterface_Internal(DM);
331: /* Applications may use this function */
332: PETSC_EXTERN PetscErrorCode DMPlexCreateNumbering_Plex(DM, PetscInt, PetscInt, PetscInt, PetscInt *, PetscSF, IS *);
334: PETSC_INTERN PetscErrorCode DMPlexCreateCellNumbering_Internal(DM, PetscBool, IS *);
335: PETSC_INTERN PetscErrorCode DMPlexCreateVertexNumbering_Internal(DM, PetscBool, IS *);
336: PETSC_INTERN PetscErrorCode DMPlexRefine_Internal(DM, DMLabel, DM *);
337: PETSC_INTERN PetscErrorCode DMPlexCoarsen_Internal(DM, DMLabel, DM *);
338: PETSC_INTERN PetscErrorCode DMCreateMatrix_Plex(DM, Mat*);
340: PETSC_INTERN PetscErrorCode DMPlexGetOverlap_Plex(DM, PetscInt *);
342: /* invert dihedral symmetry: return a^-1,
343: * using the representation described in
344: * DMPlexGetConeOrientation() */
345: PETSC_STATIC_INLINE PetscInt DihedralInvert(PetscInt N, PetscInt a)
346: {
347: return (a <= 0) ? a : (N - a);
348: }
350: /* invert dihedral symmetry: return b * a,
351: * using the representation described in
352: * DMPlexGetConeOrientation() */
353: PETSC_STATIC_INLINE PetscInt DihedralCompose(PetscInt N, PetscInt a, PetscInt b)
354: {
355: if (!N) return 0;
356: return (a >= 0) ?
357: ((b >= 0) ? ((a + b) % N) : -(((a - b - 1) % N) + 1)) :
358: ((b >= 0) ? -(((N - b - a - 1) % N) + 1) : ((N + b - a) % N));
359: }
361: /* swap dihedral symmetries: return b * a^-1,
362: * using the representation described in
363: * DMPlexGetConeOrientation() */
364: PETSC_STATIC_INLINE PetscInt DihedralSwap(PetscInt N, PetscInt a, PetscInt b)
365: {
366: return DihedralCompose(N,DihedralInvert(N,a),b);
367: }
369: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Internal(DM, IS , PetscReal, Vec, Vec, PetscReal, Vec, void *);
370: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Internal(DM, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
371: PETSC_EXTERN PetscErrorCode DMPlexReconstructGradients_Internal(DM, PetscFV, PetscInt, PetscInt, Vec, Vec, Vec, Vec);
373: /* Matvec with A in row-major storage, x and y can be aliased */
374: PETSC_STATIC_INLINE void DMPlex_Mult2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
375: {
376: PetscScalar z[2];
377: z[0] = x[0]; z[1] = x[ldx];
378: y[0] = A[0]*z[0] + A[1]*z[1];
379: y[ldx] = A[2]*z[0] + A[3]*z[1];
380: (void)PetscLogFlops(6.0);
381: }
382: PETSC_STATIC_INLINE void DMPlex_Mult3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
383: {
384: PetscScalar z[3];
385: z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
386: y[0] = A[0]*z[0] + A[1]*z[1] + A[2]*z[2];
387: y[ldx] = A[3]*z[0] + A[4]*z[1] + A[5]*z[2];
388: y[ldx*2] = A[6]*z[0] + A[7]*z[1] + A[8]*z[2];
389: (void)PetscLogFlops(15.0);
390: }
391: PETSC_STATIC_INLINE void DMPlex_MultTranspose2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
392: {
393: PetscScalar z[2];
394: z[0] = x[0]; z[1] = x[ldx];
395: y[0] = A[0]*z[0] + A[2]*z[1];
396: y[ldx] = A[1]*z[0] + A[3]*z[1];
397: (void)PetscLogFlops(6.0);
398: }
399: PETSC_STATIC_INLINE void DMPlex_MultTranspose3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
400: {
401: PetscScalar z[3];
402: z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
403: y[0] = A[0]*z[0] + A[3]*z[1] + A[6]*z[2];
404: y[ldx] = A[1]*z[0] + A[4]*z[1] + A[7]*z[2];
405: y[ldx*2] = A[2]*z[0] + A[5]*z[1] + A[8]*z[2];
406: (void)PetscLogFlops(15.0);
407: }
408: PETSC_STATIC_INLINE void DMPlex_Mult2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
409: {
410: PetscScalar z[2];
411: z[0] = x[0]; z[1] = x[ldx];
412: y[0] = A[0]*z[0] + A[1]*z[1];
413: y[ldx] = A[2]*z[0] + A[3]*z[1];
414: (void)PetscLogFlops(6.0);
415: }
416: PETSC_STATIC_INLINE void DMPlex_Mult3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
417: {
418: PetscScalar z[3];
419: z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
420: y[0] = A[0]*z[0] + A[1]*z[1] + A[2]*z[2];
421: y[ldx] = A[3]*z[0] + A[4]*z[1] + A[5]*z[2];
422: y[ldx*2] = A[6]*z[0] + A[7]*z[1] + A[8]*z[2];
423: (void)PetscLogFlops(15.0);
424: }
425: PETSC_STATIC_INLINE void DMPlex_MultTranspose2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
426: {
427: PetscScalar z[2];
428: z[0] = x[0]; z[1] = x[ldx];
429: y[0] = A[0]*z[0] + A[2]*z[1];
430: y[ldx] = A[1]*z[0] + A[3]*z[1];
431: (void)PetscLogFlops(6.0);
432: }
433: PETSC_STATIC_INLINE void DMPlex_MultTranspose3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
434: {
435: PetscScalar z[3];
436: z[0] = x[0]; z[1] = x[ldx]; z[2] = x[ldx*2];
437: y[0] = A[0]*z[0] + A[3]*z[1] + A[6]*z[2];
438: y[ldx] = A[1]*z[0] + A[4]*z[1] + A[7]*z[2];
439: y[ldx*2] = A[2]*z[0] + A[5]*z[1] + A[8]*z[2];
440: (void)PetscLogFlops(15.0);
441: }
443: PETSC_STATIC_INLINE void DMPlex_MatMult2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
444: {
445: PetscInt j;
446: for (j = 0; j < n; ++j) {
447: PetscScalar z[2];
448: z[0] = B[0+j]; z[1] = B[1*ldb+j];
449: DMPlex_Mult2D_Internal(A, 1, z, z);
450: C[0+j] = z[0]; C[1*ldb+j] = z[1];
451: }
452: (void)PetscLogFlops(8.0*n);
453: }
454: PETSC_STATIC_INLINE void DMPlex_MatMult3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
455: {
456: PetscInt j;
457: for (j = 0; j < n; ++j) {
458: PetscScalar z[3];
459: z[0] = B[0+j]; z[1] = B[1*ldb+j]; z[2] = B[2*ldb+j];
460: DMPlex_Mult3D_Internal(A, 1, z, z);
461: C[0+j] = z[0]; C[1*ldb+j] = z[1]; C[2*ldb+j] = z[2];
462: }
463: (void)PetscLogFlops(8.0*n);
464: }
465: PETSC_STATIC_INLINE void DMPlex_MatMultTranspose2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
466: {
467: PetscInt j;
468: for (j = 0; j < n; ++j) {
469: PetscScalar z[2];
470: z[0] = B[0+j]; z[1] = B[1*ldb+j];
471: DMPlex_MultTranspose2D_Internal(A, 1, z, z);
472: C[0+j] = z[0]; C[1*ldb+j] = z[1];
473: }
474: (void)PetscLogFlops(8.0*n);
475: }
476: PETSC_STATIC_INLINE void DMPlex_MatMultTranspose3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
477: {
478: PetscInt j;
479: for (j = 0; j < n; ++j) {
480: PetscScalar z[3];
481: z[0] = B[0+j]; z[1] = B[1*ldb+j]; z[2] = B[2*ldb+j];
482: DMPlex_MultTranspose3D_Internal(A, 1, z, z);
483: C[0+j] = z[0]; C[1*ldb+j] = z[1]; C[2*ldb+j] = z[2];
484: }
485: (void)PetscLogFlops(8.0*n);
486: }
488: PETSC_STATIC_INLINE void DMPlex_MatMultLeft2D_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_MultTranspose2D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
493: }
494: (void)PetscLogFlops(8.0*m);
495: }
496: PETSC_STATIC_INLINE void DMPlex_MatMultLeft3D_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_MultTranspose3D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
501: }
502: (void)PetscLogFlops(8.0*m);
503: }
504: PETSC_STATIC_INLINE void DMPlex_MatMultTransposeLeft2D_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_Mult2D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
509: }
510: (void)PetscLogFlops(8.0*m);
511: }
512: PETSC_STATIC_INLINE void DMPlex_MatMultTransposeLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
513: {
514: PetscInt j;
515: for (j = 0; j < m; ++j) {
516: DMPlex_Mult3D_Internal(A, 1, &B[j*ldb], &C[j*ldb]);
517: }
518: (void)PetscLogFlops(8.0*m);
519: }
521: PETSC_STATIC_INLINE void DMPlex_Transpose2D_Internal(PetscScalar A[])
522: {
523: PetscScalar tmp;
524: tmp = A[1]; A[1] = A[2]; A[2] = tmp;
525: }
526: PETSC_STATIC_INLINE void DMPlex_Transpose3D_Internal(PetscScalar A[])
527: {
528: PetscScalar tmp;
529: tmp = A[1]; A[1] = A[3]; A[3] = tmp;
530: tmp = A[2]; A[2] = A[6]; A[6] = tmp;
531: tmp = A[5]; A[5] = A[7]; A[7] = tmp;
532: }
534: PETSC_STATIC_INLINE void DMPlex_Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
535: {
536: const PetscReal invDet = 1.0/detJ;
538: invJ[0] = invDet*J[3];
539: invJ[1] = -invDet*J[1];
540: invJ[2] = -invDet*J[2];
541: invJ[3] = invDet*J[0];
542: (void)PetscLogFlops(5.0);
543: }
545: PETSC_STATIC_INLINE void DMPlex_Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
546: {
547: const PetscReal invDet = 1.0/detJ;
549: invJ[0*3+0] = invDet*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]);
550: invJ[0*3+1] = invDet*(J[0*3+2]*J[2*3+1] - J[0*3+1]*J[2*3+2]);
551: invJ[0*3+2] = invDet*(J[0*3+1]*J[1*3+2] - J[0*3+2]*J[1*3+1]);
552: invJ[1*3+0] = invDet*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]);
553: invJ[1*3+1] = invDet*(J[0*3+0]*J[2*3+2] - J[0*3+2]*J[2*3+0]);
554: invJ[1*3+2] = invDet*(J[0*3+2]*J[1*3+0] - J[0*3+0]*J[1*3+2]);
555: invJ[2*3+0] = invDet*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]);
556: invJ[2*3+1] = invDet*(J[0*3+1]*J[2*3+0] - J[0*3+0]*J[2*3+1]);
557: invJ[2*3+2] = invDet*(J[0*3+0]*J[1*3+1] - J[0*3+1]*J[1*3+0]);
558: (void)PetscLogFlops(37.0);
559: }
561: PETSC_STATIC_INLINE void DMPlex_Det2D_Internal(PetscReal *detJ, const PetscReal J[])
562: {
563: *detJ = J[0]*J[3] - J[1]*J[2];
564: (void)PetscLogFlops(3.0);
565: }
567: PETSC_STATIC_INLINE void DMPlex_Det3D_Internal(PetscReal *detJ, const PetscReal J[])
568: {
569: *detJ = (J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
570: J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
571: J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
572: (void)PetscLogFlops(12.0);
573: }
575: PETSC_STATIC_INLINE void DMPlex_Det2D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
576: {
577: *detJ = PetscRealPart(J[0])*PetscRealPart(J[3]) - PetscRealPart(J[1])*PetscRealPart(J[2]);
578: (void)PetscLogFlops(3.0);
579: }
581: PETSC_STATIC_INLINE void DMPlex_Det3D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
582: {
583: *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])) +
584: 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])) +
585: 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])));
586: (void)PetscLogFlops(12.0);
587: }
589: 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];}
591: 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;}
593: 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;}
595: 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);}
597: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Translate_Private(PetscInt ornt, PetscInt *start, PetscBool *reverse)
598: {
600: *reverse = (ornt < 0) ? PETSC_TRUE : PETSC_FALSE;
601: *start = *reverse ? -(ornt+1) : ornt;
602: return(0);
603: }
605: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Combine_Private(PetscInt coneSize, PetscInt origStart, PetscBool origReverse, PetscInt rotateStart, PetscBool rotateReverse, PetscInt *newStart, PetscBool *newReverse)
606: {
608: *newReverse = (origReverse == rotateReverse) ? PETSC_FALSE : PETSC_TRUE;
609: *newStart = rotateReverse ? (coneSize + rotateStart - origStart) : (coneSize + origStart - rotateStart);
610: *newStart %= coneSize;
611: return(0);
612: }
614: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_TranslateBack_Private(PetscInt coneSize, PetscInt start, PetscBool reverse, PetscInt *ornt)
615: {
617: if (coneSize < 3) {
618: /* edges just get flipped if start == 1 regardless direction */
619: *ornt = start ? -2 : 0;
620: } else {
621: *ornt = reverse ? -(start+1) : start;
622: }
623: return(0);
624: }
626: PETSC_STATIC_INLINE PetscErrorCode DMPlexFixFaceOrientations_Permute_Private(PetscInt n, const PetscInt arr[], PetscInt start, PetscBool reverse, PetscInt newarr[])
627: {
628: PetscInt i;
631: if (reverse) {for (i=0; i<n; i++) newarr[i] = arr[(n+start-i)%n];}
632: else {for (i=0; i<n; i++) newarr[i] = arr[(start+i)%n];}
633: return(0);
634: }
636: PETSC_INTERN PetscErrorCode DMPlexGetPointDualSpaceFEM(DM,PetscInt,PetscInt,PetscDualSpace *);
637: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPoint_Internal(PetscSection,PetscBool,PetscInt,PetscInt,PetscInt *,PetscBool,const PetscInt[],const PetscInt[],PetscInt[]);
638: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPointFields_Internal(PetscSection,PetscBool,PetscInt,PetscInt,PetscInt[],PetscBool,const PetscInt***,PetscInt,const PetscInt[],PetscInt[]);
639: PETSC_INTERN PetscErrorCode DMPlexGetCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);
640: PETSC_INTERN PetscErrorCode DMPlexRestoreCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);
642: PETSC_EXTERN PetscErrorCode DMSNESGetFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
643: PETSC_EXTERN PetscErrorCode DMSNESRestoreFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
644: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Patch_Internal(DM, PetscSection, IS, PetscReal, Vec, Vec, Vec, void *);
645: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Patch_Internal(DM, PetscSection, PetscSection, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
646: PETSC_INTERN PetscErrorCode DMCreateSubDomainDM_Plex(DM,DMLabel,PetscInt,IS*,DM*);
647: PETSC_INTERN PetscErrorCode DMPlexBasisTransformPoint_Internal(DM, DM, Vec, PetscInt, PetscBool[], PetscBool, PetscScalar *);
648: PETSC_EXTERN PetscErrorCode DMPlexBasisTransformPointTensor_Internal(DM, DM, Vec, PetscInt, PetscBool, PetscInt, PetscScalar *);
649: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApplyReal_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscReal *, PetscReal *, void *);
650: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApply_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscScalar *, PetscScalar *, void *);
651: PETSC_INTERN PetscErrorCode DMCreateNeumannOverlap_Plex(DM, IS*, Mat*, PetscErrorCode (**)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void*), void **);
653: #endif /* _PLEXIMPL_H */