Actual source code: dmpleximpl.h
1: #pragma once
3: #include <petscmat.h>
4: #include <petscdmplex.h>
5: #include <petscdmplextransform.h>
6: #include <petscbt.h>
7: #include <petscsf.h>
8: #include <petsc/private/dmimpl.h>
10: #if defined(PETSC_HAVE_EXODUSII)
11: #include <exodusII.h>
12: #endif
14: PETSC_EXTERN PetscBool Plexcite;
15: PETSC_EXTERN const char PlexCitation[];
17: PETSC_EXTERN PetscLogEvent DMPLEX_Interpolate;
18: PETSC_EXTERN PetscLogEvent DMPLEX_Partition;
19: PETSC_EXTERN PetscLogEvent DMPLEX_PartSelf;
20: PETSC_EXTERN PetscLogEvent DMPLEX_PartLabelInvert;
21: PETSC_EXTERN PetscLogEvent DMPLEX_PartLabelCreateSF;
22: PETSC_EXTERN PetscLogEvent DMPLEX_PartStratSF;
23: PETSC_EXTERN PetscLogEvent DMPLEX_CreatePointSF;
24: PETSC_EXTERN PetscLogEvent DMPLEX_Distribute;
25: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeCones;
26: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeLabels;
27: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeSF;
28: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeOverlap;
29: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeField;
30: PETSC_EXTERN PetscLogEvent DMPLEX_DistributeData;
31: PETSC_EXTERN PetscLogEvent DMPLEX_Migrate;
32: PETSC_EXTERN PetscLogEvent DMPLEX_InterpolateSF;
33: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalToNaturalBegin;
34: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalToNaturalEnd;
35: PETSC_EXTERN PetscLogEvent DMPLEX_NaturalToGlobalBegin;
36: PETSC_EXTERN PetscLogEvent DMPLEX_NaturalToGlobalEnd;
37: PETSC_EXTERN PetscLogEvent DMPLEX_Stratify;
38: PETSC_EXTERN PetscLogEvent DMPLEX_Symmetrize;
39: PETSC_EXTERN PetscLogEvent DMPLEX_Preallocate;
40: PETSC_EXTERN PetscLogEvent DMPLEX_ResidualFEM;
41: PETSC_EXTERN PetscLogEvent DMPLEX_JacobianFEM;
42: PETSC_EXTERN PetscLogEvent DMPLEX_InterpolatorFEM;
43: PETSC_EXTERN PetscLogEvent DMPLEX_InjectorFEM;
44: PETSC_EXTERN PetscLogEvent DMPLEX_IntegralFEM;
45: PETSC_EXTERN PetscLogEvent DMPLEX_CreateGmsh;
46: PETSC_EXTERN PetscLogEvent DMPLEX_RebalanceSharedPoints;
47: PETSC_EXTERN PetscLogEvent DMPLEX_CreateFromFile;
48: PETSC_EXTERN PetscLogEvent DMPLEX_CreateFromOptions;
49: PETSC_EXTERN PetscLogEvent DMPLEX_BuildFromCellList;
50: PETSC_EXTERN PetscLogEvent DMPLEX_BuildCoordinatesFromCellList;
51: PETSC_EXTERN PetscLogEvent DMPLEX_LocatePoints;
52: PETSC_EXTERN PetscLogEvent DMPLEX_TopologyView;
53: PETSC_EXTERN PetscLogEvent DMPLEX_DistributionView;
54: PETSC_EXTERN PetscLogEvent DMPLEX_LabelsView;
55: PETSC_EXTERN PetscLogEvent DMPLEX_CoordinatesView;
56: PETSC_EXTERN PetscLogEvent DMPLEX_SectionView;
57: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalVectorView;
58: PETSC_EXTERN PetscLogEvent DMPLEX_LocalVectorView;
59: PETSC_EXTERN PetscLogEvent DMPLEX_TopologyLoad;
60: PETSC_EXTERN PetscLogEvent DMPLEX_DistributionLoad;
61: PETSC_EXTERN PetscLogEvent DMPLEX_LabelsLoad;
62: PETSC_EXTERN PetscLogEvent DMPLEX_CoordinatesLoad;
63: PETSC_EXTERN PetscLogEvent DMPLEX_SectionLoad;
64: PETSC_EXTERN PetscLogEvent DMPLEX_GlobalVectorLoad;
65: PETSC_EXTERN PetscLogEvent DMPLEX_LocalVectorLoad;
66: PETSC_EXTERN PetscLogEvent DMPLEX_MetricEnforceSPD;
67: PETSC_EXTERN PetscLogEvent DMPLEX_MetricNormalize;
68: PETSC_EXTERN PetscLogEvent DMPLEX_MetricAverage;
69: PETSC_EXTERN PetscLogEvent DMPLEX_MetricIntersection;
70: PETSC_EXTERN PetscLogEvent DMPLEX_Generate;
71: PETSC_EXTERN PetscLogEvent DMPLEX_Transform;
72: PETSC_EXTERN PetscLogEvent DMPLEX_GetLocalOffsets;
74: PETSC_EXTERN PetscLogEvent DMPLEX_RebalBuildGraph;
75: PETSC_EXTERN PetscLogEvent DMPLEX_RebalRewriteSF;
76: PETSC_EXTERN PetscLogEvent DMPLEX_RebalGatherGraph;
77: PETSC_EXTERN PetscLogEvent DMPLEX_RebalPartition;
78: PETSC_EXTERN PetscLogEvent DMPLEX_RebalScatterPart;
79: PETSC_EXTERN PetscLogEvent DMPLEX_Uninterpolate;
81: /* Utility struct to store the contents of a Fluent file in memory */
82: typedef struct {
83: int index; /* Type of section */
84: unsigned int zoneID;
85: unsigned int first;
86: unsigned int last;
87: int type;
88: int nd; /* Either ND or element-type */
89: void *data;
90: } FluentSection;
92: struct _PetscGridHash {
93: PetscInt dim;
94: PetscReal lower[3]; /* The lower-left corner */
95: PetscReal upper[3]; /* The upper-right corner */
96: PetscReal extent[3]; /* The box size */
97: PetscReal h[3]; /* The subbox size */
98: PetscInt n[3]; /* The number of subboxes */
99: PetscSection cellSection; /* Offsets for cells in each subbox*/
100: IS cells; /* List of cells in each subbox */
101: DMLabel cellsSparse; /* Sparse storage for cell map */
102: };
104: typedef struct {
105: PetscBool isotropic; /* Is the metric isotropic? */
106: PetscBool uniform; /* Is the metric uniform? */
107: PetscBool restrictAnisotropyFirst; /* Should anisotropy or normalization come first? */
108: PetscBool noInsert; /* Should node insertion/deletion be turned off? */
109: PetscBool noSwap; /* Should facet swapping be turned off? */
110: PetscBool noMove; /* Should node movement be turned off? */
111: PetscBool noSurf; /* Should surface modification be turned off? */
112: PetscReal h_min, h_max; /* Minimum/maximum tolerated metric magnitudes */
113: PetscReal a_max; /* Maximum tolerated anisotropy */
114: PetscReal targetComplexity; /* Target metric complexity */
115: PetscReal p; /* Degree for L-p normalization methods */
116: PetscReal gradationFactor; /* Maximum tolerated length ratio for adjacent edges */
117: PetscReal hausdorffNumber; /* Max. distance between piecewise linear representation of boundary and reconstructed ideal boundary */
118: PetscInt numIter; /* Number of ParMmg mesh adaptation iterations */
119: PetscInt verbosity; /* Level of verbosity for remesher (-1 = no output, 10 = maximum) */
120: } DMPlexMetricCtx;
122: /* Point Numbering in Plex:
124: Points are numbered contiguously by stratum. Strate are organized as follows:
126: First Stratum: Cells [height 0]
127: Second Stratum: Vertices [depth 0]
128: Third Stratum: Faces [height 1]
129: Fourth Stratum: Edges [depth 1]
131: We do this so that the numbering of a cell-vertex mesh does not change after interpolation. Within a given stratum,
132: we allow additional segregation of by cell type.
133: */
134: typedef struct {
135: PetscInt refct;
137: PetscSection coneSection; /* Layout of cones (inedges for DAG) */
138: PetscInt *cones; /* Cone for each point */
139: PetscInt *coneOrientations; /* Orientation of each cone point, means cone traversal should start on point 'o', and if negative start on -(o+1) and go in reverse */
140: PetscSection supportSection; /* Layout of cones (inedges for DAG) */
141: PetscInt *supports; /* Cone for each point */
143: struct { // DMPolytopeType is an enum (usually size 4), but this needs frequent access
144: uint8_t value_as_uint8; // in a struct to guard for stronger typing
145: } *cellTypes;
147: /* Transformation */
148: DMPlexTransform tr; /* Type of transform used to define an ephemeral mesh */
149: char *transformType; /* Type of transform for uniform cell refinement */
150: PetscBool refinementUniform; /* Flag for uniform cell refinement */
151: PetscReal refinementLimit; /* Maximum volume for refined cell */
152: PetscErrorCode (*refinementFunc)(const PetscReal[], PetscReal *); /* Function giving the maximum volume for refined cell */
154: /* Interpolation */
155: DMPlexInterpolatedFlag interpolated;
156: DMPlexInterpolatedFlag interpolatedCollective;
158: DMPlexReorderDefaultFlag reorderDefault; /* Reorder the DM by default */
160: /* Distribution */
161: PetscBool distDefault; /* Distribute the DM by default */
162: PetscInt overlap; /* Overlap of the partitions as passed to DMPlexDistribute() or DMPlexDistributeOverlap() */
163: PetscInt numOvLabels; /* The number of labels used for candidate overlap points */
164: DMLabel ovLabels[16]; /* Labels used for candidate overlap points */
165: PetscInt ovValues[16]; /* Label values used for candidate overlap points */
166: PetscInt numOvExLabels; /* The number of labels used for exclusion */
167: DMLabel ovExLabels[16]; /* Labels used to exclude points from the overlap */
168: PetscInt ovExValues[16]; /* Label values used to exclude points from the overlap */
169: char *distributionName; /* Name of the specific parallel distribution of the DM */
171: /* Hierarchy */
172: PetscBool regularRefinement; /* This flag signals that we are a regular refinement of coarseMesh */
174: /* Generation */
175: char *tetgenOpts;
176: char *triangleOpts;
177: PetscPartitioner partitioner;
178: PetscBool partitionBalance; /* Evenly divide partition overlap when distributing */
179: PetscBool remeshBd;
181: /* Submesh */
182: DMLabel subpointMap; /* Label each original mesh point in the submesh with its depth, subpoint are the implicit numbering */
183: IS subpointIS; /* IS holding point number in the enclosing mesh of every point in the submesh chart */
184: PetscObjectState subpointState; /* The state of subpointMap when the subpointIS was last created */
186: /* Labels and numbering */
187: PetscObjectState depthState; /* State of depth label, so that we can determine if a user changes it */
188: PetscObjectState celltypeState; /* State of celltype label, so that we can determine if a user changes it */
189: IS globalVertexNumbers;
190: IS globalCellNumbers;
192: /* Constraints */
193: PetscSection anchorSection; /* maps constrained points to anchor points */
194: IS anchorIS; /* anchors indexed by the above section */
195: PetscErrorCode (*createanchors)(DM); /* automatically compute anchors (probably from tree constraints) */
196: PetscErrorCode (*computeanchormatrix)(DM, PetscSection, PetscSection, Mat);
198: /* Tree: automatically construct constraints for hierarchically non-conforming meshes */
199: PetscSection parentSection; /* dof == 1 if point has parent */
200: PetscInt *parents; /* point to parent */
201: PetscInt *childIDs; /* point to child ID */
202: PetscSection childSection; /* inverse of parent section */
203: PetscInt *children; /* point to children */
204: DM referenceTree; /* reference tree to which child ID's refer */
205: PetscErrorCode (*getchildsymmetry)(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *);
207: /* MATIS support */
208: PetscSection subdomainSection;
210: /* Adjacency */
211: PetscBool useAnchors; /* Replace constrained points with their anchors in adjacency lists */
212: PetscErrorCode (*useradjacency)(DM, PetscInt, PetscInt *, PetscInt[], void *); /* User callback for adjacency */
213: void *useradjacencyctx; /* User context for callback */
215: // Periodicity
216: struct {
217: // Specified by the user
218: PetscScalar transform[4][4]; // geometric transform
219: PetscSF face_sf; // root(donor faces) <-- leaf(local faces)
220: // Created eagerly (depends on points)
221: PetscSF composed_sf; // root(non-periodic global points) <-- leaf(local points)
222: IS periodic_points;
223: } periodic;
225: /* Projection */
226: PetscInt maxProjectionHeight; /* maximum height of cells used in DMPlexProject functions */
227: PetscInt activePoint; /* current active point in iteration */
229: /* Output */
230: PetscInt vtkCellHeight; /* The height of cells for output, default is 0 */
231: PetscReal scale[NUM_PETSC_UNITS]; /* The scale for each SI unit */
233: /* Geometry */
234: PetscBool ignoreModel; /* Ignore the geometry model during refinement */
235: PetscReal minradius; /* Minimum distance from cell centroid to face */
236: PetscBool useHashLocation; /* Use grid hashing for point location */
237: PetscGridHash lbox; /* Local box for searching */
238: void (*coordFunc)(PetscInt, PetscInt, PetscInt, /* Function used to remap newly introduced vertices */
239: 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[]);
241: /* Neighbors */
242: PetscMPIInt *neighbors;
244: /* Metric */
245: DMPlexMetricCtx *metricCtx;
247: /* FEM */
248: PetscBool useCeed; /* This should convert to a registration system when there are more FEM backends */
249: PetscBool useMatClPerm; /* Use the closure permutation when assembling matrices */
251: /* Debugging */
252: PetscBool printSetValues;
253: PetscInt printFEM;
254: PetscInt printL2;
255: PetscInt printLocate;
256: PetscReal printTol;
257: } DM_Plex;
259: PETSC_INTERN PetscErrorCode DMPlexCopy_Internal(DM, PetscBool, PetscBool, DM);
260: PETSC_INTERN PetscErrorCode DMPlexReplace_Internal(DM, DM *);
262: PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll_VTU(DM, PetscViewer);
263: PETSC_EXTERN PetscErrorCode VecView_Plex_Local(Vec, PetscViewer);
264: PETSC_EXTERN PetscErrorCode VecView_Plex_Native(Vec, PetscViewer);
265: PETSC_EXTERN PetscErrorCode VecView_Plex(Vec, PetscViewer);
266: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Local(Vec, PetscViewer);
267: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Native(Vec, PetscViewer);
268: PETSC_EXTERN PetscErrorCode VecLoad_Plex(Vec, PetscViewer);
269: PETSC_INTERN PetscErrorCode DMPlexGetFieldType_Internal(DM, PetscSection, PetscInt, PetscInt *, PetscInt *, PetscViewerVTKFieldType *);
270: PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM, PetscViewer);
271: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject, PetscViewer);
272: #if defined(PETSC_HAVE_HDF5)
273: PETSC_EXTERN PetscErrorCode VecView_Plex_Local_HDF5(Vec, PetscViewer);
274: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5(Vec, PetscViewer);
275: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5(Vec, PetscViewer);
276: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5_Native(Vec, PetscViewer);
277: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5_Native(Vec, PetscViewer);
278: PETSC_EXTERN PetscErrorCode DMPlexView_HDF5(DM, PetscViewer);
279: PETSC_EXTERN PetscErrorCode DMPlexLoad_HDF5(DM, PetscViewer);
280: PETSC_INTERN PetscErrorCode DMPlexTopologyView_HDF5_Internal(DM, IS, PetscViewer);
281: PETSC_INTERN PetscErrorCode DMPlexCoordinatesView_HDF5_Internal(DM, PetscViewer);
282: PETSC_INTERN PetscErrorCode DMPlexLabelsView_HDF5_Internal(DM, IS, PetscViewer);
283: PETSC_INTERN PetscErrorCode DMPlexSectionView_HDF5_Internal(DM, PetscViewer, DM);
284: PETSC_INTERN PetscErrorCode DMPlexGlobalVectorView_HDF5_Internal(DM, PetscViewer, DM, Vec);
285: PETSC_INTERN PetscErrorCode DMPlexLocalVectorView_HDF5_Internal(DM, PetscViewer, DM, Vec);
286: PETSC_INTERN PetscErrorCode DMPlexTopologyLoad_HDF5_Internal(DM, PetscViewer, PetscSF *);
287: PETSC_INTERN PetscErrorCode DMPlexCoordinatesLoad_HDF5_Internal(DM, PetscViewer, PetscSF);
288: PETSC_INTERN PetscErrorCode DMPlexLabelsLoad_HDF5_Internal(DM, PetscViewer, PetscSF);
289: PETSC_INTERN PetscErrorCode DMPlexSectionLoad_HDF5_Internal(DM, PetscViewer, DM, PetscSF, PetscSF *, PetscSF *);
290: PETSC_INTERN PetscErrorCode DMPlexVecLoad_HDF5_Internal(DM, PetscViewer, DM, PetscSF, Vec);
291: PETSC_INTERN PetscErrorCode DMPlexView_HDF5_Internal(DM, PetscViewer);
292: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Internal(DM, PetscViewer);
293: PETSC_INTERN PetscErrorCode DMPlexLoad_HDF5_Xdmf_Internal(DM, PetscViewer);
294: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Internal(Vec, PetscViewer);
295: PETSC_INTERN PetscErrorCode VecView_Plex_HDF5_Native_Internal(Vec, PetscViewer);
296: PETSC_INTERN PetscErrorCode VecView_Plex_Local_HDF5_Internal(Vec, PetscViewer);
297: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Internal(Vec, PetscViewer);
298: PETSC_INTERN PetscErrorCode VecLoad_Plex_HDF5_Native_Internal(Vec, PetscViewer);
300: struct _n_DMPlexStorageVersion {
301: int major, minor, subminor;
302: };
303: typedef struct _n_DMPlexStorageVersion *DMPlexStorageVersion;
305: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetDMPlexStorageVersionReading(PetscViewer, DMPlexStorageVersion *);
306: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetDMPlexStorageVersionWriting(PetscViewer, DMPlexStorageVersion *);
307: #endif
308: PETSC_EXTERN PetscErrorCode VecView_Plex_Local_CGNS(Vec, PetscViewer);
310: PETSC_INTERN PetscErrorCode DMPlexVecGetClosure_Internal(DM, PetscSection, PetscBool, Vec, PetscInt, PetscInt *, PetscScalar *[]);
311: PETSC_INTERN PetscErrorCode DMPlexVecGetClosureAtDepth_Internal(DM, PetscSection, Vec, PetscInt, PetscInt, PetscInt *, PetscScalar *[]);
312: PETSC_INTERN PetscErrorCode DMPlexClosurePoints_Private(DM, PetscInt, const PetscInt[], IS *);
313: PETSC_INTERN PetscErrorCode DMSetFromOptions_NonRefinement_Plex(DM, PetscOptionItems *);
314: PETSC_INTERN PetscErrorCode DMSetFromOptions_Overlap_Plex(DM, PetscOptionItems *, PetscInt *);
315: PETSC_INTERN PetscErrorCode DMCoarsen_Plex(DM, MPI_Comm, DM *);
316: PETSC_INTERN PetscErrorCode DMCoarsenHierarchy_Plex(DM, PetscInt, DM[]);
317: PETSC_INTERN PetscErrorCode DMRefine_Plex(DM, MPI_Comm, DM *);
318: PETSC_INTERN PetscErrorCode DMRefineHierarchy_Plex(DM, PetscInt, DM[]);
319: PETSC_INTERN PetscErrorCode DMAdaptLabel_Plex(DM, Vec, DMLabel, DMLabel, DM *);
320: PETSC_INTERN PetscErrorCode DMExtrude_Plex(DM, PetscInt, DM *);
321: PETSC_INTERN PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
322: PETSC_INTERN PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
323: PETSC_INTERN PetscErrorCode DMProjectFunctionLocal_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
324: 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);
325: 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);
326: 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);
327: PETSC_INTERN PetscErrorCode DMProjectBdFieldLabelLocal_Plex(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Vec, void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode, Vec);
328: PETSC_INTERN PetscErrorCode DMComputeL2Diff_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
329: PETSC_INTERN PetscErrorCode DMComputeL2GradientDiff_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal[], PetscReal *);
330: PETSC_INTERN PetscErrorCode DMComputeL2FieldDiff_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
331: PETSC_INTERN PetscErrorCode DMLocatePoints_Plex(DM, Vec, DMPointLocationType, PetscSF);
333: #if defined(PETSC_HAVE_EXODUSII)
334: PETSC_INTERN PetscErrorCode DMView_PlexExodusII(DM, PetscViewer);
335: PETSC_INTERN PetscErrorCode VecView_PlexExodusII_Internal(Vec, PetscViewer);
336: PETSC_INTERN PetscErrorCode VecLoad_PlexExodusII_Internal(Vec, PetscViewer);
337: #endif
338: PETSC_INTERN PetscErrorCode DMView_PlexCGNS(DM, PetscViewer);
339: PETSC_INTERN PetscErrorCode DMPlexCreateCGNSFromFile_Internal(MPI_Comm, const char *, PetscBool, DM *);
340: PETSC_INTERN PetscErrorCode DMPlexCreateCGNS_Internal(MPI_Comm, PetscInt, PetscBool, DM *);
341: PETSC_INTERN PetscErrorCode DMPlexVTKGetCellType_Internal(DM, PetscInt, PetscInt, PetscInt *);
342: PETSC_INTERN PetscErrorCode DMPlexGetAdjacency_Internal(DM, PetscInt, PetscBool, PetscBool, PetscBool, PetscInt *, PetscInt *[]);
343: PETSC_INTERN PetscErrorCode DMPlexGetRawFaces_Internal(DM, DMPolytopeType, const PetscInt[], PetscInt *, const DMPolytopeType *[], const PetscInt *[], const PetscInt *[]);
344: PETSC_INTERN PetscErrorCode DMPlexRestoreRawFaces_Internal(DM, DMPolytopeType, const PetscInt[], PetscInt *, const DMPolytopeType *[], const PetscInt *[], const PetscInt *[]);
345: PETSC_INTERN PetscErrorCode DMPlexComputeCellType_Internal(DM, PetscInt, PetscInt, DMPolytopeType *);
346: PETSC_INTERN PetscErrorCode DMPlexVecSetFieldClosure_Internal(DM, PetscSection, Vec, PetscBool[], PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscScalar[], InsertMode);
347: PETSC_INTERN PetscErrorCode DMPlexProjectConstraints_Internal(DM, Vec, Vec);
348: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_SetTree(DM, PetscSection, PetscInt[], PetscInt[]);
349: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_Union(DM, DM, const char *, DM *);
350: PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorTree(DM, DM, PetscSF, PetscInt *, Mat);
351: PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorTree(DM, DM, PetscSF, PetscInt *, Mat);
352: PETSC_EXTERN PetscErrorCode DMPlexAnchorsModifyMat(DM, PetscSection, PetscInt, PetscInt, const PetscInt[], const PetscInt ***, const PetscScalar[], PetscInt *, PetscInt *, PetscInt *[], PetscScalar *[], PetscInt[], PetscBool);
353: PETSC_INTERN PetscErrorCode DMPlexLocatePoint_Internal(DM, PetscInt, const PetscScalar[], PetscInt, PetscInt *);
354: /* these two are PETSC_EXTERN just because of src/dm/impls/plex/tests/ex18.c */
355: PETSC_EXTERN PetscErrorCode DMPlexOrientInterface_Internal(DM);
357: /* Applications may use this function */
358: PETSC_EXTERN PetscErrorCode DMPlexCreateNumbering_Plex(DM, PetscInt, PetscInt, PetscInt, PetscInt *, PetscSF, IS *);
360: PETSC_INTERN PetscErrorCode DMPlexInterpolateInPlace_Internal(DM);
361: PETSC_INTERN PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Internal(DM, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool);
362: PETSC_INTERN PetscErrorCode DMPlexMigrateIsoperiodicFaceSF_Internal(DM, DM, PetscSF);
363: PETSC_INTERN PetscErrorCode DMPlexCreateCellNumbering_Internal(DM, PetscBool, IS *);
364: PETSC_INTERN PetscErrorCode DMPlexCreateVertexNumbering_Internal(DM, PetscBool, IS *);
365: PETSC_INTERN PetscErrorCode DMPlexRefine_Internal(DM, Vec, DMLabel, DMLabel, DM *);
366: PETSC_INTERN PetscErrorCode DMPlexCoarsen_Internal(DM, Vec, DMLabel, DMLabel, DM *);
367: PETSC_INTERN PetscErrorCode DMCreateMatrix_Plex(DM, Mat *);
369: PETSC_INTERN PetscErrorCode DMPlexGetOverlap_Plex(DM, PetscInt *);
370: PETSC_INTERN PetscErrorCode DMPlexSetOverlap_Plex(DM, DM, PetscInt);
371: PETSC_INTERN PetscErrorCode DMPlexDistributeGetDefault_Plex(DM, PetscBool *);
372: PETSC_INTERN PetscErrorCode DMPlexDistributeSetDefault_Plex(DM, PetscBool);
373: PETSC_INTERN PetscErrorCode DMPlexReorderGetDefault_Plex(DM, DMPlexReorderDefaultFlag *);
374: PETSC_INTERN PetscErrorCode DMPlexReorderSetDefault_Plex(DM, DMPlexReorderDefaultFlag);
375: PETSC_INTERN PetscErrorCode DMPlexGetUseCeed_Plex(DM, PetscBool *);
376: PETSC_INTERN PetscErrorCode DMPlexSetUseCeed_Plex(DM, PetscBool);
378: #if 1
379: static inline PetscInt DihedralInvert(PetscInt N, PetscInt a)
380: {
381: return (a <= 0) ? a : (N - a);
382: }
384: static inline PetscInt DihedralCompose(PetscInt N, PetscInt a, PetscInt b)
385: {
386: if (!N) return 0;
387: return (a >= 0) ? ((b >= 0) ? ((a + b) % N) : -(((a - b - 1) % N) + 1)) : ((b >= 0) ? -(((N - b - a - 1) % N) + 1) : ((N + b - a) % N));
388: }
390: static inline PetscInt DihedralSwap(PetscInt N, PetscInt a, PetscInt b)
391: {
392: return DihedralCompose(N, DihedralInvert(N, a), b);
393: }
394: #else
395: /* TODO
396: This is a reimplementation of the tensor dihedral symmetries using the new orientations.
397: These should be turned on when we convert to new-style orientations in p4est.
398: */
399: /* invert dihedral symmetry: return a^-1,
400: * using the representation described in
401: * DMPlexGetConeOrientation() */
402: static inline PetscInt DihedralInvert(PetscInt N, PetscInt a)
403: {
404: switch (N) {
405: case 0:
406: return 0;
407: case 2:
408: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_SEGMENT, 0, a);
409: case 4:
410: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_QUADRILATERAL, 0, a);
411: case 8:
412: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_HEXAHEDRON, 0, a);
413: default:
414: SETERRABORT(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid celltype for DihedralInvert()");
415: }
416: return 0;
417: }
419: /* compose dihedral symmetry: return b * a,
420: * using the representation described in
421: * DMPlexGetConeOrientation() */
422: static inline PetscInt DihedralCompose(PetscInt N, PetscInt a, PetscInt b)
423: {
424: switch (N) {
425: case 0:
426: return 0;
427: case 2:
428: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_SEGMENT, b, a);
429: case 4:
430: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_QUADRILATERAL, b, a);
431: case 8:
432: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_HEXAHEDRON, b, a);
433: default:
434: SETERRABORT(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid celltype for DihedralCompose()");
435: }
436: return 0;
437: }
439: /* swap dihedral symmetries: return b * a^-1,
440: * using the representation described in
441: * DMPlexGetConeOrientation() */
442: static inline PetscInt DihedralSwap(PetscInt N, PetscInt a, PetscInt b)
443: {
444: switch (N) {
445: case 0:
446: return 0;
447: case 2:
448: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_SEGMENT, b, a);
449: case 4:
450: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_QUADRILATERAL, b, a);
451: case 8:
452: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_HEXAHEDRON, b, a);
453: default:
454: SETERRABORT(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid celltype for DihedralCompose()");
455: }
456: return 0;
457: }
458: #endif
459: PETSC_EXTERN PetscInt DMPolytopeConvertNewOrientation_Internal(DMPolytopeType, PetscInt);
460: PETSC_EXTERN PetscInt DMPolytopeConvertOldOrientation_Internal(DMPolytopeType, PetscInt);
461: PETSC_EXTERN PetscErrorCode DMPlexConvertOldOrientations_Internal(DM);
463: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Internal(DM, PetscFormKey, IS, PetscReal, Vec, Vec, PetscReal, Vec, void *);
464: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Hybrid_Internal(DM, PetscFormKey[], IS, PetscReal, Vec, Vec, PetscReal, Vec, void *);
465: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Internal(DM, PetscFormKey, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
466: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Hybrid_Internal(DM, PetscFormKey[], IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
467: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Action_Internal(DM, PetscFormKey, IS, PetscReal, PetscReal, Vec, Vec, Vec, Vec, void *);
468: PETSC_EXTERN PetscErrorCode DMPlexReconstructGradients_Internal(DM, PetscFV, PetscInt, PetscInt, Vec, Vec, Vec, Vec);
470: /* Matvec with A in row-major storage, x and y can be aliased */
471: static inline void DMPlex_Mult2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
472: {
473: const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
474: y[0 * ldx] = A[0] * z[0] + A[1] * z[1];
475: y[1 * ldx] = A[2] * z[0] + A[3] * z[1];
476: (void)PetscLogFlops(6.0);
477: }
478: static inline void DMPlex_Mult3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
479: {
480: const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
481: y[0 * ldx] = A[0] * z[0] + A[1] * z[1] + A[2] * z[2];
482: y[1 * ldx] = A[3] * z[0] + A[4] * z[1] + A[5] * z[2];
483: y[2 * ldx] = A[6] * z[0] + A[7] * z[1] + A[8] * z[2];
484: (void)PetscLogFlops(15.0);
485: }
486: static inline void DMPlex_MultTranspose2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
487: {
488: const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
489: y[0 * ldx] = A[0] * z[0] + A[2] * z[1];
490: y[1 * ldx] = A[1] * z[0] + A[3] * z[1];
491: (void)PetscLogFlops(6.0);
492: }
493: static inline void DMPlex_MultTranspose3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
494: {
495: const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
496: y[0 * ldx] = PetscConj(A[0]) * z[0] + PetscConj(A[3]) * z[1] + PetscConj(A[6]) * z[2];
497: y[1 * ldx] = PetscConj(A[1]) * z[0] + PetscConj(A[4]) * z[1] + PetscConj(A[7]) * z[2];
498: y[2 * ldx] = PetscConj(A[2]) * z[0] + PetscConj(A[5]) * z[1] + PetscConj(A[8]) * z[2];
499: (void)PetscLogFlops(15.0);
500: }
501: static inline void DMPlex_Mult2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
502: {
503: const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
504: y[0 * ldx] = A[0] * z[0] + A[1] * z[1];
505: y[1 * ldx] = A[2] * z[0] + A[3] * z[1];
506: (void)PetscLogFlops(6.0);
507: }
508: static inline void DMPlex_Mult3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
509: {
510: const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
511: y[0 * ldx] = A[0] * z[0] + A[1] * z[1] + A[2] * z[2];
512: y[1 * ldx] = A[3] * z[0] + A[4] * z[1] + A[5] * z[2];
513: y[2 * ldx] = A[6] * z[0] + A[7] * z[1] + A[8] * z[2];
514: (void)PetscLogFlops(15.0);
515: }
516: static inline void DMPlex_MultAdd2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
517: {
518: const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
519: y[0 * ldx] += A[0] * z[0] + A[1] * z[1];
520: y[1 * ldx] += A[2] * z[0] + A[3] * z[1];
521: (void)PetscLogFlops(6.0);
522: }
523: static inline void DMPlex_MultAdd3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
524: {
525: const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
526: y[0 * ldx] += A[0] * z[0] + A[1] * z[1] + A[2] * z[2];
527: y[1 * ldx] += A[3] * z[0] + A[4] * z[1] + A[5] * z[2];
528: y[2 * ldx] += A[6] * z[0] + A[7] * z[1] + A[8] * z[2];
529: (void)PetscLogFlops(15.0);
530: }
531: /*
532: A: packed, row-major m x n array
533: x: length m array
534: y: length n arra
535: ldx: the stride in x and y
537: A[i,j] = A[i * n + j]
538: A^T[j,i] = A[i,j]
539: */
540: static inline void DMPlex_MultTransposeReal_Internal(const PetscReal A[], PetscInt m, PetscInt n, PetscInt ldx, const PetscScalar x[], PetscScalar y[])
541: {
542: PetscScalar z[3];
543: PetscInt i, j;
544: for (i = 0; i < m; ++i) z[i] = x[i * ldx];
545: for (j = 0; j < n; ++j) {
546: const PetscInt l = j * ldx;
547: y[l] = 0;
548: for (i = 0; i < m; ++i) y[l] += A[i * n + j] * z[i];
549: }
550: (void)PetscLogFlops(2 * m * n);
551: }
552: static inline void DMPlex_MultTranspose2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
553: {
554: const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
555: y[0 * ldx] = A[0] * z[0] + A[2] * z[1];
556: y[1 * ldx] = A[1] * z[0] + A[3] * z[1];
557: (void)PetscLogFlops(6.0);
558: }
559: static inline void DMPlex_MultTranspose3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
560: {
561: const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
562: y[0 * ldx] = A[0] * z[0] + A[3] * z[1] + A[6] * z[2];
563: y[1 * ldx] = A[1] * z[0] + A[4] * z[1] + A[7] * z[2];
564: y[2 * ldx] = A[2] * z[0] + A[5] * z[1] + A[8] * z[2];
565: (void)PetscLogFlops(15.0);
566: }
568: static inline void DMPlex_MatMult2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
569: {
570: #define PLEX_DIM__ 2
571: PetscScalar z[PLEX_DIM__];
572: for (PetscInt j = 0; j < n; ++j) {
573: for (int d = 0; d < PLEX_DIM__; ++d) z[d] = B[d * ldb + j];
574: DMPlex_Mult2D_Internal(A, 1, z, z);
575: for (int d = 0; d < PLEX_DIM__; ++d) C[d * ldb + j] = z[d];
576: }
577: (void)PetscLogFlops(8.0 * n);
578: #undef PLEX_DIM__
579: }
580: static inline void DMPlex_MatMult3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
581: {
582: #define PLEX_DIM__ 3
583: PetscScalar z[PLEX_DIM__];
584: for (PetscInt j = 0; j < n; ++j) {
585: for (int d = 0; d < PLEX_DIM__; ++d) z[d] = B[d * ldb + j];
586: DMPlex_Mult3D_Internal(A, 1, z, z);
587: for (int d = 0; d < PLEX_DIM__; ++d) C[d * ldb + j] = z[d];
588: }
589: (void)PetscLogFlops(8.0 * n);
590: #undef PLEX_DIM__
591: }
592: static inline void DMPlex_MatMultTranspose2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
593: {
594: #define PLEX_DIM__ 2
595: PetscScalar z[PLEX_DIM__];
596: for (PetscInt j = 0; j < n; ++j) {
597: for (int d = 0; d < PLEX_DIM__; ++d) z[d] = B[d * ldb + j];
598: DMPlex_MultTranspose2D_Internal(A, 1, z, z);
599: for (int d = 0; d < PLEX_DIM__; ++d) C[d * ldb + j] = z[d];
600: }
601: (void)PetscLogFlops(8.0 * n);
602: #undef PLEX_DIM__
603: }
604: static inline void DMPlex_MatMultTranspose3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
605: {
606: #define PLEX_DIM__ 3
607: PetscScalar z[PLEX_DIM__];
608: for (PetscInt j = 0; j < n; ++j) {
609: for (int d = 0; d < PLEX_DIM__; ++d) z[d] = B[d * ldb + j];
610: DMPlex_MultTranspose3D_Internal(A, 1, z, z);
611: for (int d = 0; d < PLEX_DIM__; ++d) C[d * ldb + j] = z[d];
612: }
613: (void)PetscLogFlops(8.0 * n);
614: #undef PLEX_DIM__
615: }
617: static inline void DMPlex_MatMultLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
618: {
619: for (PetscInt j = 0; j < m; ++j) DMPlex_MultTranspose2D_Internal(A, 1, &B[j * ldb], &C[j * ldb]);
620: (void)PetscLogFlops(8.0 * m);
621: }
622: static inline void DMPlex_MatMultLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
623: {
624: for (PetscInt j = 0; j < m; ++j) DMPlex_MultTranspose3D_Internal(A, 1, &B[j * ldb], &C[j * ldb]);
625: (void)PetscLogFlops(8.0 * m);
626: }
627: static inline void DMPlex_MatMultTransposeLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
628: {
629: for (PetscInt j = 0; j < m; ++j) DMPlex_Mult2D_Internal(A, 1, &B[j * ldb], &C[j * ldb]);
630: (void)PetscLogFlops(8.0 * m);
631: }
632: static inline void DMPlex_MatMultTransposeLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
633: {
634: for (PetscInt j = 0; j < m; ++j) DMPlex_Mult3D_Internal(A, 1, &B[j * ldb], &C[j * ldb]);
635: (void)PetscLogFlops(8.0 * m);
636: }
638: static inline void DMPlex_PTAP2DReal_Internal(const PetscReal P[], const PetscScalar A[], PetscScalar C[])
639: {
640: PetscScalar out[4];
641: PetscInt i, j, k, l;
642: for (i = 0; i < 2; ++i) {
643: for (j = 0; j < 2; ++j) {
644: out[i * 2 + j] = 0.;
645: for (k = 0; k < 2; ++k) {
646: for (l = 0; l < 2; ++l) out[i * 2 + j] += P[k * 2 + i] * A[k * 2 + l] * P[l * 2 + j];
647: }
648: }
649: }
650: for (i = 0; i < 2 * 2; ++i) C[i] = out[i];
651: (void)PetscLogFlops(48.0);
652: }
653: static inline void DMPlex_PTAP3DReal_Internal(const PetscReal P[], const PetscScalar A[], PetscScalar C[])
654: {
655: PetscScalar out[9];
656: PetscInt i, j, k, l;
657: for (i = 0; i < 3; ++i) {
658: for (j = 0; j < 3; ++j) {
659: out[i * 3 + j] = 0.;
660: for (k = 0; k < 3; ++k) {
661: for (l = 0; l < 3; ++l) out[i * 3 + j] += P[k * 3 + i] * A[k * 3 + l] * P[l * 3 + j];
662: }
663: }
664: }
665: for (i = 0; i < 2 * 2; ++i) C[i] = out[i];
666: (void)PetscLogFlops(243.0);
667: }
668: /* TODO Fix for aliasing of A and C */
669: static inline void DMPlex_PTAPReal_Internal(const PetscReal P[], PetscInt m, PetscInt n, const PetscScalar A[], PetscScalar C[])
670: {
671: PetscInt i, j, k, l;
672: for (i = 0; i < n; ++i) {
673: for (j = 0; j < n; ++j) {
674: C[i * n + j] = 0.;
675: for (k = 0; k < m; ++k) {
676: for (l = 0; l < m; ++l) C[i * n + j] += P[k * n + i] * A[k * m + l] * P[l * n + j];
677: }
678: }
679: }
680: (void)PetscLogFlops(243.0);
681: }
683: static inline void DMPlex_Transpose2D_Internal(PetscScalar A[])
684: {
685: PetscScalar tmp;
686: tmp = A[1];
687: A[1] = A[2];
688: A[2] = tmp;
689: }
690: static inline void DMPlex_Transpose3D_Internal(PetscScalar A[])
691: {
692: PetscScalar tmp;
693: tmp = A[1];
694: A[1] = A[3];
695: A[3] = tmp;
696: tmp = A[2];
697: A[2] = A[6];
698: A[6] = tmp;
699: tmp = A[5];
700: A[5] = A[7];
701: A[7] = tmp;
702: }
704: static inline void DMPlex_Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
705: {
706: // Allow zero volume cells
707: const PetscReal invDet = detJ == 0 ? 1.0 : 1.0 / detJ;
709: invJ[0] = invDet * J[3];
710: invJ[1] = -invDet * J[1];
711: invJ[2] = -invDet * J[2];
712: invJ[3] = invDet * J[0];
713: (void)PetscLogFlops(5.0);
714: }
716: static inline void DMPlex_Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
717: {
718: // Allow zero volume cells
719: const PetscReal invDet = detJ == 0 ? 1.0 : 1.0 / detJ;
721: invJ[0 * 3 + 0] = invDet * (J[1 * 3 + 1] * J[2 * 3 + 2] - J[1 * 3 + 2] * J[2 * 3 + 1]);
722: invJ[0 * 3 + 1] = invDet * (J[0 * 3 + 2] * J[2 * 3 + 1] - J[0 * 3 + 1] * J[2 * 3 + 2]);
723: invJ[0 * 3 + 2] = invDet * (J[0 * 3 + 1] * J[1 * 3 + 2] - J[0 * 3 + 2] * J[1 * 3 + 1]);
724: invJ[1 * 3 + 0] = invDet * (J[1 * 3 + 2] * J[2 * 3 + 0] - J[1 * 3 + 0] * J[2 * 3 + 2]);
725: invJ[1 * 3 + 1] = invDet * (J[0 * 3 + 0] * J[2 * 3 + 2] - J[0 * 3 + 2] * J[2 * 3 + 0]);
726: invJ[1 * 3 + 2] = invDet * (J[0 * 3 + 2] * J[1 * 3 + 0] - J[0 * 3 + 0] * J[1 * 3 + 2]);
727: invJ[2 * 3 + 0] = invDet * (J[1 * 3 + 0] * J[2 * 3 + 1] - J[1 * 3 + 1] * J[2 * 3 + 0]);
728: invJ[2 * 3 + 1] = invDet * (J[0 * 3 + 1] * J[2 * 3 + 0] - J[0 * 3 + 0] * J[2 * 3 + 1]);
729: invJ[2 * 3 + 2] = invDet * (J[0 * 3 + 0] * J[1 * 3 + 1] - J[0 * 3 + 1] * J[1 * 3 + 0]);
730: (void)PetscLogFlops(37.0);
731: }
733: static inline void DMPlex_Det2D_Internal(PetscReal *detJ, const PetscReal J[])
734: {
735: *detJ = J[0] * J[3] - J[1] * J[2];
736: (void)PetscLogFlops(3.0);
737: }
739: static inline void DMPlex_Det3D_Internal(PetscReal *detJ, const PetscReal J[])
740: {
741: *detJ = (J[0 * 3 + 0] * (J[1 * 3 + 1] * J[2 * 3 + 2] - J[1 * 3 + 2] * J[2 * 3 + 1]) + J[0 * 3 + 1] * (J[1 * 3 + 2] * J[2 * 3 + 0] - J[1 * 3 + 0] * J[2 * 3 + 2]) + J[0 * 3 + 2] * (J[1 * 3 + 0] * J[2 * 3 + 1] - J[1 * 3 + 1] * J[2 * 3 + 0]));
742: (void)PetscLogFlops(12.0);
743: }
745: static inline void DMPlex_Det2D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
746: {
747: *detJ = PetscRealPart(J[0]) * PetscRealPart(J[3]) - PetscRealPart(J[1]) * PetscRealPart(J[2]);
748: (void)PetscLogFlops(3.0);
749: }
751: static inline void DMPlex_Det3D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
752: {
753: *detJ = (PetscRealPart(J[0 * 3 + 0]) * (PetscRealPart(J[1 * 3 + 1]) * PetscRealPart(J[2 * 3 + 2]) - PetscRealPart(J[1 * 3 + 2]) * PetscRealPart(J[2 * 3 + 1])) + PetscRealPart(J[0 * 3 + 1]) * (PetscRealPart(J[1 * 3 + 2]) * PetscRealPart(J[2 * 3 + 0]) - PetscRealPart(J[1 * 3 + 0]) * PetscRealPart(J[2 * 3 + 2])) + PetscRealPart(J[0 * 3 + 2]) * (PetscRealPart(J[1 * 3 + 0]) * PetscRealPart(J[2 * 3 + 1]) - PetscRealPart(J[1 * 3 + 1]) * PetscRealPart(J[2 * 3 + 0])));
754: (void)PetscLogFlops(12.0);
755: }
757: static inline void DMPlex_WaxpyD_Internal(PetscInt dim, PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w)
758: {
759: PetscInt d;
760: for (d = 0; d < dim; ++d) w[d] = a * x[d] + y[d];
761: }
763: static inline PetscReal DMPlex_DotD_Internal(PetscInt dim, const PetscScalar *x, const PetscReal *y)
764: {
765: PetscReal sum = 0.0;
766: PetscInt d;
767: for (d = 0; d < dim; ++d) sum += PetscRealPart(x[d]) * y[d];
768: return sum;
769: }
771: static inline PetscReal DMPlex_DotRealD_Internal(PetscInt dim, const PetscReal *x, const PetscReal *y)
772: {
773: PetscReal sum = 0.0;
774: PetscInt d;
775: for (d = 0; d < dim; ++d) sum += x[d] * y[d];
776: return sum;
777: }
779: static inline PetscReal DMPlex_NormD_Internal(PetscInt dim, const PetscReal *x)
780: {
781: PetscReal sum = 0.0;
782: PetscInt d;
783: for (d = 0; d < dim; ++d) sum += x[d] * x[d];
784: return PetscSqrtReal(sum);
785: }
787: static inline PetscReal DMPlex_DistD_Internal(PetscInt dim, const PetscScalar *x, const PetscScalar *y)
788: {
789: PetscReal sum = 0.0;
790: PetscInt d;
791: for (d = 0; d < dim; ++d) sum += PetscRealPart(PetscConj(x[d] - y[d]) * (x[d] - y[d]));
792: return PetscSqrtReal(sum);
793: }
795: PETSC_INTERN PetscErrorCode DMPlexGetPointDualSpaceFEM(DM, PetscInt, PetscInt, PetscDualSpace *);
796: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPoint_Internal(PetscSection, PetscBool, PetscInt, PetscInt, PetscInt *, PetscBool, const PetscInt[], const PetscInt[], PetscInt[]);
797: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPointFields_Internal(PetscSection, PetscBool, PetscInt, PetscInt, PetscInt[], PetscBool, const PetscInt ***, PetscInt, const PetscInt[], PetscInt[]);
798: PETSC_INTERN PetscErrorCode DMPlexGetTransitiveClosure_Internal(DM, PetscInt, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
799: PETSC_INTERN PetscErrorCode DMPlexVecGetOrientedClosure_Internal(DM, PetscSection, PetscBool, Vec, PetscInt, PetscInt, PetscInt *, PetscScalar *[]);
800: PETSC_INTERN PetscErrorCode DMPlexMatSetClosure_Internal(DM, PetscSection, PetscSection, PetscBool, Mat, PetscInt, const PetscScalar[], InsertMode);
802: PETSC_EXTERN PetscErrorCode DMPlexGetAllCells_Internal(DM, IS *);
803: PETSC_EXTERN PetscErrorCode DMSNESGetFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
804: PETSC_EXTERN PetscErrorCode DMSNESRestoreFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
805: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Patch_Internal(DM, PetscSection, IS, PetscReal, Vec, Vec, Vec, void *);
806: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Patch_Internal(DM, PetscSection, PetscSection, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
807: PETSC_INTERN PetscErrorCode DMCreateSubDomainDM_Plex(DM, DMLabel, PetscInt, IS *, DM *);
808: PETSC_INTERN PetscErrorCode DMPlexBasisTransformPoint_Internal(DM, DM, Vec, PetscInt, PetscBool[], PetscBool, PetscScalar *);
809: PETSC_EXTERN PetscErrorCode DMPlexBasisTransformPointTensor_Internal(DM, DM, Vec, PetscInt, PetscBool, PetscInt, PetscScalar *);
810: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApplyReal_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscReal *, PetscReal *, void *);
811: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApply_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscScalar *, PetscScalar *, void *);
812: PETSC_INTERN PetscErrorCode DMCreateNeumannOverlap_Plex(DM, IS *, Mat *, PetscErrorCode (**)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void *), void **);
814: PETSC_INTERN PetscErrorCode DMPeriodicCoordinateSetUp_Internal(DM);
816: /* Functions in the vtable */
817: PETSC_INTERN PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
818: PETSC_INTERN PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat);
819: PETSC_INTERN PetscErrorCode DMCreateMassMatrix_Plex(DM dmCoarse, DM dmFine, Mat *mat);
820: PETSC_INTERN PetscErrorCode DMCreateMassMatrixLumped_Plex(DM dmCoarse, Vec *lm);
821: PETSC_INTERN PetscErrorCode DMCreateLocalSection_Plex(DM dm);
822: PETSC_INTERN PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm);
823: PETSC_INTERN PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J);
824: PETSC_INTERN PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
825: PETSC_INTERN PetscErrorCode DMCreateCoordinateField_Plex(DM dm, DMField *field);
826: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
827: PETSC_INTERN PetscErrorCode DMSetUp_Plex(DM dm);
828: PETSC_INTERN PetscErrorCode DMDestroy_Plex(DM dm);
829: PETSC_INTERN PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
830: PETSC_INTERN PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
831: PETSC_INTERN PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm);
832: PETSC_INTERN PetscErrorCode DMCreateSuperDM_Plex(DM dms[], PetscInt len, IS **is, DM *superdm);