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: /* Ordering */
159: DMReorderDefaultFlag reorderDefault; /* Reorder the DM by default */
161: /* Distribution */
162: PetscBool distDefault; /* Distribute the DM by default */
163: PetscInt overlap; /* Overlap of the partitions as passed to DMPlexDistribute() or DMPlexDistributeOverlap() */
164: PetscInt numOvLabels; /* The number of labels used for candidate overlap points */
165: DMLabel ovLabels[16]; /* Labels used for candidate overlap points */
166: PetscInt ovValues[16]; /* Label values used for candidate overlap points */
167: PetscInt numOvExLabels; /* The number of labels used for exclusion */
168: DMLabel ovExLabels[16]; /* Labels used to exclude points from the overlap */
169: PetscInt ovExValues[16]; /* Label values used to exclude points from the overlap */
170: char *distributionName; /* Name of the specific parallel distribution of the DM */
172: /* Hierarchy */
173: PetscBool regularRefinement; /* This flag signals that we are a regular refinement of coarseMesh */
175: /* Generation */
176: char *tetgenOpts;
177: char *triangleOpts;
178: PetscPartitioner partitioner;
179: PetscBool partitionBalance; /* Evenly divide partition overlap when distributing */
180: PetscBool remeshBd;
182: /* Submesh */
183: DMLabel subpointMap; /* Label each original mesh point in the submesh with its depth, subpoint are the implicit numbering */
184: IS subpointIS; /* IS holding point number in the enclosing mesh of every point in the submesh chart */
185: PetscObjectState subpointState; /* The state of subpointMap when the subpointIS was last created */
187: /* Labels and numbering */
188: PetscObjectState depthState; /* State of depth label, so that we can determine if a user changes it */
189: PetscObjectState celltypeState; /* State of celltype label, so that we can determine if a user changes it */
190: IS globalVertexNumbers;
191: IS globalCellNumbers;
193: /* Constraints */
194: PetscSection anchorSection; /* maps constrained points to anchor points */
195: IS anchorIS; /* anchors indexed by the above section */
196: PetscErrorCode (*createanchors)(DM); /* automatically compute anchors (probably from tree constraints) */
197: PetscErrorCode (*computeanchormatrix)(DM, PetscSection, PetscSection, Mat);
199: /* Tree: automatically construct constraints for hierarchically non-conforming meshes */
200: PetscSection parentSection; /* dof == 1 if point has parent */
201: PetscInt *parents; /* point to parent */
202: PetscInt *childIDs; /* point to child ID */
203: PetscSection childSection; /* inverse of parent section */
204: PetscInt *children; /* point to children */
205: DM referenceTree; /* reference tree to which child ID's refer */
206: PetscErrorCode (*getchildsymmetry)(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *);
208: /* MATIS support */
209: PetscSection subdomainSection;
211: /* Adjacency */
212: PetscBool useAnchors; /* Replace constrained points with their anchors in adjacency lists */
213: PetscErrorCode (*useradjacency)(DM, PetscInt, PetscInt *, PetscInt[], void *); /* User callback for adjacency */
214: void *useradjacencyctx; /* User context for callback */
216: // Periodicity
217: struct {
218: // Specified by the user
219: PetscInt num_face_sfs; // number of face_sfs
220: PetscSF *face_sfs; // root(donor faces) <-- leaf(local faces)
221: PetscScalar (*transform)[4][4]; // geometric transform
222: // Created eagerly (depends on points)
223: PetscSF composed_sf; // root(non-periodic global points) <-- leaf(local points)
224: IS *periodic_points;
225: } periodic;
227: /* Projection */
228: PetscInt maxProjectionHeight; /* maximum height of cells used in DMPlexProject functions */
229: PetscInt activePoint; /* current active point in iteration */
231: /* Output */
232: PetscInt vtkCellHeight; /* The height of cells for output, default is 0 */
233: PetscReal scale[NUM_PETSC_UNITS]; /* The scale for each SI unit */
235: /* Geometry */
236: PetscBool ignoreModel; /* Ignore the geometry model during refinement */
237: PetscReal minradius; /* Minimum distance from cell centroid to face */
238: PetscBool useHashLocation; /* Use grid hashing for point location */
239: PetscGridHash lbox; /* Local box for searching */
240: void (*coordFunc)(PetscInt, PetscInt, PetscInt, /* Function used to remap newly introduced vertices */
241: 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[]);
243: /* Neighbors */
244: PetscMPIInt *neighbors;
246: /* Metric */
247: DMPlexMetricCtx *metricCtx;
249: /* FEM */
250: PetscBool useCeed; /* This should convert to a registration system when there are more FEM backends */
251: PetscBool useMatClPerm; /* Use the closure permutation when assembling matrices */
253: /* Debugging */
254: PetscBool printSetValues;
255: PetscInt printFEM;
256: PetscInt printFVM;
257: PetscInt printL2;
258: PetscInt printLocate;
259: PetscReal printTol;
260: } DM_Plex;
262: PETSC_INTERN PetscErrorCode DMPlexCopy_Internal(DM, PetscBool, PetscBool, DM);
263: PETSC_INTERN PetscErrorCode DMPlexReplace_Internal(DM, DM *);
265: PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll_VTU(DM, PetscViewer);
266: PETSC_EXTERN PetscErrorCode VecView_Plex_Local(Vec, PetscViewer);
267: PETSC_EXTERN PetscErrorCode VecView_Plex_Native(Vec, PetscViewer);
268: PETSC_EXTERN PetscErrorCode VecView_Plex(Vec, PetscViewer);
269: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Local(Vec, PetscViewer);
270: PETSC_EXTERN PetscErrorCode VecLoad_Plex_Native(Vec, PetscViewer);
271: PETSC_EXTERN PetscErrorCode VecLoad_Plex(Vec, PetscViewer);
272: PETSC_INTERN PetscErrorCode DMPlexGetFieldTypes_Internal(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscInt **, PetscViewerVTKFieldType **);
273: PETSC_INTERN PetscErrorCode DMPlexRestoreFieldTypes_Internal(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscInt **, PetscViewerVTKFieldType **);
274: PETSC_INTERN PetscErrorCode DMPlexGetFieldType_Internal(DM, PetscSection, PetscInt, PetscInt *, PetscInt *, PetscViewerVTKFieldType *);
275: PETSC_INTERN PetscErrorCode DMPlexView_GLVis(DM, PetscViewer);
276: PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_Plex(PetscObject, PetscViewer);
277: #if defined(PETSC_HAVE_HDF5)
278: PETSC_EXTERN PetscErrorCode VecView_Plex_Local_HDF5(Vec, PetscViewer);
279: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5(Vec, PetscViewer);
280: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5(Vec, PetscViewer);
281: PETSC_EXTERN PetscErrorCode VecView_Plex_HDF5_Native(Vec, PetscViewer);
282: PETSC_EXTERN PetscErrorCode VecLoad_Plex_HDF5_Native(Vec, PetscViewer);
283: PETSC_EXTERN PetscErrorCode DMPlexView_HDF5(DM, PetscViewer);
284: PETSC_EXTERN PetscErrorCode DMPlexLoad_HDF5(DM, PetscViewer);
285: PETSC_INTERN PetscErrorCode DMPlexTopologyView_HDF5_Internal(DM, IS, PetscViewer);
286: PETSC_INTERN PetscErrorCode DMPlexCoordinatesView_HDF5_Internal(DM, PetscViewer);
287: PETSC_INTERN PetscErrorCode DMPlexLabelsView_HDF5_Internal(DM, IS, PetscViewer);
288: PETSC_INTERN PetscErrorCode DMPlexSectionView_HDF5_Internal(DM, PetscViewer, DM);
289: PETSC_INTERN PetscErrorCode DMPlexGlobalVectorView_HDF5_Internal(DM, PetscViewer, DM, Vec);
290: PETSC_INTERN PetscErrorCode DMPlexLocalVectorView_HDF5_Internal(DM, PetscViewer, DM, Vec);
291: PETSC_INTERN PetscErrorCode DMPlexTopologyLoad_HDF5_Internal(DM, PetscViewer, PetscSF *);
292: PETSC_INTERN PetscErrorCode DMPlexCoordinatesLoad_HDF5_Internal(DM, PetscViewer, PetscSF);
293: PETSC_INTERN PetscErrorCode DMPlexLabelsLoad_HDF5_Internal(DM, PetscViewer, PetscSF);
294: PETSC_INTERN PetscErrorCode DMPlexSectionLoad_HDF5_Internal(DM, PetscViewer, DM, PetscSF, PetscSF *, PetscSF *);
295: PETSC_INTERN PetscErrorCode DMPlexVecLoad_HDF5_Internal(DM, PetscViewer, DM, PetscSF, Vec);
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);
305: struct _n_DMPlexStorageVersion {
306: int major, minor, subminor;
307: };
308: typedef struct _n_DMPlexStorageVersion *DMPlexStorageVersion;
310: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetDMPlexStorageVersionReading(PetscViewer, DMPlexStorageVersion *);
311: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetDMPlexStorageVersionWriting(PetscViewer, DMPlexStorageVersion *);
312: #endif
313: PETSC_EXTERN PetscErrorCode VecView_Plex_Local_CGNS(Vec, PetscViewer);
315: PETSC_INTERN PetscErrorCode DMPlexVecGetClosure_Internal(DM, PetscSection, PetscBool, Vec, PetscInt, PetscInt *, PetscScalar *[]);
316: PETSC_INTERN PetscErrorCode DMPlexVecGetClosureAtDepth_Internal(DM, PetscSection, Vec, PetscInt, PetscInt, PetscInt *, PetscScalar *[]);
317: PETSC_INTERN PetscErrorCode DMPlexClosurePoints_Private(DM, PetscInt, const PetscInt[], IS *);
318: PETSC_INTERN PetscErrorCode DMSetFromOptions_NonRefinement_Plex(DM, PetscOptionItems *);
319: PETSC_INTERN PetscErrorCode DMSetFromOptions_Overlap_Plex(DM, PetscOptionItems *, PetscInt *);
320: PETSC_INTERN PetscErrorCode DMCoarsen_Plex(DM, MPI_Comm, DM *);
321: PETSC_INTERN PetscErrorCode DMCoarsenHierarchy_Plex(DM, PetscInt, DM[]);
322: PETSC_INTERN PetscErrorCode DMRefine_Plex(DM, MPI_Comm, DM *);
323: PETSC_INTERN PetscErrorCode DMRefineHierarchy_Plex(DM, PetscInt, DM[]);
324: PETSC_INTERN PetscErrorCode DMAdaptLabel_Plex(DM, Vec, DMLabel, DMLabel, DM *);
325: PETSC_INTERN PetscErrorCode DMExtrude_Plex(DM, PetscInt, DM *);
326: PETSC_INTERN PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
327: PETSC_INTERN PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues_Plex(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
328: PETSC_INTERN PetscErrorCode DMProjectFunctionLocal_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
329: 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);
330: 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);
331: 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);
332: 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);
333: PETSC_INTERN PetscErrorCode DMComputeL2Diff_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
334: PETSC_INTERN PetscErrorCode DMComputeL2GradientDiff_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal[], PetscReal *);
335: PETSC_INTERN PetscErrorCode DMComputeL2FieldDiff_Plex(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
336: PETSC_INTERN PetscErrorCode DMLocatePoints_Plex(DM, Vec, DMPointLocationType, PetscSF);
338: #if defined(PETSC_HAVE_EXODUSII)
339: PETSC_INTERN PetscErrorCode DMView_PlexExodusII(DM, PetscViewer);
340: PETSC_INTERN PetscErrorCode VecView_PlexExodusII_Internal(Vec, PetscViewer);
341: PETSC_INTERN PetscErrorCode VecLoad_PlexExodusII_Internal(Vec, PetscViewer);
342: #endif
343: PETSC_INTERN PetscErrorCode DMView_PlexCGNS(DM, PetscViewer);
344: PETSC_INTERN PetscErrorCode DMPlexCreateCGNSFromFile_Internal(MPI_Comm, const char *, PetscBool, DM *);
345: PETSC_INTERN PetscErrorCode DMPlexCreateCGNS_Internal(MPI_Comm, PetscInt, PetscBool, DM *);
346: PETSC_INTERN PetscErrorCode DMPlexVTKGetCellType_Internal(DM, PetscInt, PetscInt, PetscInt *);
347: PETSC_INTERN PetscErrorCode DMPlexGetAdjacency_Internal(DM, PetscInt, PetscBool, PetscBool, PetscBool, PetscInt *, PetscInt *[]);
348: PETSC_INTERN PetscErrorCode DMPlexGetRawFaces_Internal(DM, DMPolytopeType, const PetscInt[], PetscInt *, const DMPolytopeType *[], const PetscInt *[], const PetscInt *[]);
349: PETSC_INTERN PetscErrorCode DMPlexRestoreRawFaces_Internal(DM, DMPolytopeType, const PetscInt[], PetscInt *, const DMPolytopeType *[], const PetscInt *[], const PetscInt *[]);
350: PETSC_INTERN PetscErrorCode DMPlexComputeCellType_Internal(DM, PetscInt, PetscInt, DMPolytopeType *);
351: PETSC_INTERN PetscErrorCode DMPlexVecSetFieldClosure_Internal(DM, PetscSection, Vec, PetscBool[], PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscScalar[], InsertMode);
352: PETSC_INTERN PetscErrorCode DMPlexProjectConstraints_Internal(DM, Vec, Vec);
353: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_SetTree(DM, PetscSection, PetscInt[], PetscInt[]);
354: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceTree_Union(DM, DM, const char *, DM *);
355: PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorTree(DM, DM, PetscSF, PetscInt *, Mat);
356: PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorTree(DM, DM, PetscSF, PetscInt *, Mat);
357: PETSC_INTERN PetscErrorCode DMPlexAnchorsModifyMat(DM, PetscSection, PetscInt, PetscInt, const PetscInt[], const PetscInt ***, const PetscScalar[], PetscInt *, PetscInt *, PetscInt *[], PetscScalar *[], PetscInt[], PetscBool);
358: PETSC_INTERN PetscErrorCode DMPlexAnchorsModifyMat_Internal(DM, PetscSection, PetscInt, PetscInt, const PetscInt[], const PetscInt ***, PetscInt, PetscInt, const PetscScalar[], PetscInt *, PetscInt *, PetscInt *[], PetscScalar *[], PetscInt[], PetscBool, PetscBool);
359: PETSC_INTERN PetscErrorCode DMPlexAnchorsGetSubMatModification(DM, PetscSection, PetscInt, PetscInt, const PetscInt[], const PetscInt ***, PetscInt *, PetscInt *, PetscInt *[], PetscInt[], PetscScalar *[]);
360: PETSC_INTERN PetscErrorCode DMPlexLocatePoint_Internal(DM, PetscInt, const PetscScalar[], PetscInt, PetscInt *);
361: /* this is PETSC_EXTERN just because of src/dm/impls/plex/tests/ex18.c */
362: PETSC_EXTERN PetscErrorCode DMPlexOrientInterface_Internal(DM);
363: PETSC_INTERN PetscErrorCode DMPlexOrientCells_Internal(DM, IS, IS);
365: /* Applications may use this function */
366: PETSC_EXTERN PetscErrorCode DMPlexCreateNumbering_Plex(DM, PetscInt, PetscInt, PetscInt, PetscInt *, PetscSF, IS *);
368: PETSC_INTERN PetscErrorCode DMPlexInterpolateInPlace_Internal(DM);
369: PETSC_INTERN PetscErrorCode DMPlexCreateBoxMesh_Tensor_SFC_Internal(DM, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool);
370: PETSC_INTERN PetscErrorCode DMPlexMigrateIsoperiodicFaceSF_Internal(DM, DM, PetscSF);
371: PETSC_INTERN PetscErrorCode DMPlexCreateCellNumbering_Internal(DM, PetscBool, IS *);
372: PETSC_INTERN PetscErrorCode DMPlexCreateVertexNumbering_Internal(DM, PetscBool, IS *);
373: PETSC_INTERN PetscErrorCode DMPlexRefine_Internal(DM, Vec, DMLabel, DMLabel, DM *);
374: PETSC_INTERN PetscErrorCode DMPlexCoarsen_Internal(DM, Vec, DMLabel, DMLabel, DM *);
375: PETSC_INTERN PetscErrorCode DMCreateMatrix_Plex(DM, Mat *);
377: PETSC_INTERN PetscErrorCode DMPlexGetOverlap_Plex(DM, PetscInt *);
378: PETSC_INTERN PetscErrorCode DMPlexSetOverlap_Plex(DM, DM, PetscInt);
379: PETSC_INTERN PetscErrorCode DMPlexDistributeGetDefault_Plex(DM, PetscBool *);
380: PETSC_INTERN PetscErrorCode DMPlexDistributeSetDefault_Plex(DM, PetscBool);
381: PETSC_INTERN PetscErrorCode DMPlexReorderGetDefault_Plex(DM, DMReorderDefaultFlag *);
382: PETSC_INTERN PetscErrorCode DMPlexReorderSetDefault_Plex(DM, DMReorderDefaultFlag);
383: PETSC_INTERN PetscErrorCode DMPlexGetUseCeed_Plex(DM, PetscBool *);
384: PETSC_INTERN PetscErrorCode DMPlexSetUseCeed_Plex(DM, PetscBool);
385: PETSC_INTERN PetscErrorCode DMReorderSectionGetDefault_Plex(DM, DMReorderDefaultFlag *);
386: PETSC_INTERN PetscErrorCode DMReorderSectionSetDefault_Plex(DM, DMReorderDefaultFlag);
387: PETSC_INTERN PetscErrorCode DMReorderSectionGetType_Plex(DM, MatOrderingType *);
388: PETSC_INTERN PetscErrorCode DMReorderSectionSetType_Plex(DM, MatOrderingType);
390: #if 1
391: static inline PetscInt DihedralInvert(PetscInt N, PetscInt a)
392: {
393: return (a <= 0) ? a : (N - a);
394: }
396: static inline PetscInt DihedralCompose(PetscInt N, PetscInt a, PetscInt b)
397: {
398: if (!N) return 0;
399: 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));
400: }
402: static inline PetscInt DihedralSwap(PetscInt N, PetscInt a, PetscInt b)
403: {
404: return DihedralCompose(N, DihedralInvert(N, a), b);
405: }
406: #else
407: /* TODO
408: This is a reimplementation of the tensor dihedral symmetries using the new orientations.
409: These should be turned on when we convert to new-style orientations in p4est.
410: */
411: /* invert dihedral symmetry: return a^-1,
412: * using the representation described in
413: * DMPlexGetConeOrientation() */
414: static inline PetscInt DihedralInvert(PetscInt N, PetscInt a)
415: {
416: switch (N) {
417: case 0:
418: return 0;
419: case 2:
420: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_SEGMENT, 0, a);
421: case 4:
422: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_QUADRILATERAL, 0, a);
423: case 8:
424: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_HEXAHEDRON, 0, a);
425: default:
426: SETERRABORT(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid celltype for DihedralInvert()");
427: }
428: return 0;
429: }
431: /* compose dihedral symmetry: return b * a,
432: * using the representation described in
433: * DMPlexGetConeOrientation() */
434: static inline PetscInt DihedralCompose(PetscInt N, PetscInt a, PetscInt b)
435: {
436: switch (N) {
437: case 0:
438: return 0;
439: case 2:
440: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_SEGMENT, b, a);
441: case 4:
442: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_QUADRILATERAL, b, a);
443: case 8:
444: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_HEXAHEDRON, b, a);
445: default:
446: SETERRABORT(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid celltype for DihedralCompose()");
447: }
448: return 0;
449: }
451: /* swap dihedral symmetries: return b * a^-1,
452: * using the representation described in
453: * DMPlexGetConeOrientation() */
454: static inline PetscInt DihedralSwap(PetscInt N, PetscInt a, PetscInt b)
455: {
456: switch (N) {
457: case 0:
458: return 0;
459: case 2:
460: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_SEGMENT, b, a);
461: case 4:
462: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_QUADRILATERAL, b, a);
463: case 8:
464: return DMPolytopeTypeComposeOrientationInv(DM_POLYTOPE_HEXAHEDRON, b, a);
465: default:
466: SETERRABORT(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid celltype for DihedralCompose()");
467: }
468: return 0;
469: }
470: #endif
471: PETSC_EXTERN PetscInt DMPolytopeConvertNewOrientation_Internal(DMPolytopeType, PetscInt);
472: PETSC_EXTERN PetscInt DMPolytopeConvertOldOrientation_Internal(DMPolytopeType, PetscInt);
473: PETSC_EXTERN PetscErrorCode DMPlexConvertOldOrientations_Internal(DM);
475: PETSC_EXTERN PetscErrorCode DMPlexComputeIntegral_Internal(DM, Vec, PetscInt, PetscInt, PetscScalar *, void *);
476: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Internal(DM, PetscFormKey, IS, PetscReal, Vec, Vec, PetscReal, Vec, void *);
477: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Hybrid_Internal(DM, PetscFormKey[], IS, PetscReal, Vec, Vec, PetscReal, Vec, void *);
478: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Internal(DM, PetscFormKey, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
479: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Hybrid_Internal(DM, PetscFormKey[], IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
480: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Action_Internal(DM, PetscFormKey, IS, PetscReal, PetscReal, Vec, Vec, Vec, Vec, void *);
481: PETSC_EXTERN PetscErrorCode DMPlexReconstructGradients_Internal(DM, PetscFV, PetscInt, PetscInt, Vec, Vec, Vec, Vec);
483: /* Matvec with A in row-major storage, x and y can be aliased */
484: static inline void DMPlex_Mult2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
485: {
486: const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
487: y[0 * ldx] = A[0] * z[0] + A[1] * z[1];
488: y[1 * ldx] = A[2] * z[0] + A[3] * z[1];
489: (void)PetscLogFlops(6.0);
490: }
491: static inline void DMPlex_Mult3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
492: {
493: const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
494: y[0 * ldx] = A[0] * z[0] + A[1] * z[1] + A[2] * z[2];
495: y[1 * ldx] = A[3] * z[0] + A[4] * z[1] + A[5] * z[2];
496: y[2 * ldx] = A[6] * z[0] + A[7] * z[1] + A[8] * z[2];
497: (void)PetscLogFlops(15.0);
498: }
499: static inline void DMPlex_MultTranspose2D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
500: {
501: const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
502: y[0 * ldx] = A[0] * z[0] + A[2] * z[1];
503: y[1 * ldx] = A[1] * z[0] + A[3] * z[1];
504: (void)PetscLogFlops(6.0);
505: }
506: static inline void DMPlex_MultTranspose3D_Internal(const PetscScalar A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
507: {
508: const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
509: y[0 * ldx] = PetscConj(A[0]) * z[0] + PetscConj(A[3]) * z[1] + PetscConj(A[6]) * z[2];
510: y[1 * ldx] = PetscConj(A[1]) * z[0] + PetscConj(A[4]) * z[1] + PetscConj(A[7]) * z[2];
511: y[2 * ldx] = PetscConj(A[2]) * z[0] + PetscConj(A[5]) * z[1] + PetscConj(A[8]) * z[2];
512: (void)PetscLogFlops(15.0);
513: }
514: static inline void DMPlex_Mult2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
515: {
516: const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
517: y[0 * ldx] = A[0] * z[0] + A[1] * z[1];
518: y[1 * ldx] = A[2] * z[0] + A[3] * z[1];
519: (void)PetscLogFlops(6.0);
520: }
521: static inline void DMPlex_Mult3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
522: {
523: const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
524: y[0 * ldx] = A[0] * z[0] + A[1] * z[1] + A[2] * z[2];
525: y[1 * ldx] = A[3] * z[0] + A[4] * z[1] + A[5] * z[2];
526: y[2 * ldx] = A[6] * z[0] + A[7] * z[1] + A[8] * z[2];
527: (void)PetscLogFlops(15.0);
528: }
529: static inline void DMPlex_MultAdd2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
530: {
531: const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
532: y[0 * ldx] += A[0] * z[0] + A[1] * z[1];
533: y[1 * ldx] += A[2] * z[0] + A[3] * z[1];
534: (void)PetscLogFlops(6.0);
535: }
536: static inline void DMPlex_MultAdd3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
537: {
538: const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
539: y[0 * ldx] += A[0] * z[0] + A[1] * z[1] + A[2] * z[2];
540: y[1 * ldx] += A[3] * z[0] + A[4] * z[1] + A[5] * z[2];
541: y[2 * ldx] += A[6] * z[0] + A[7] * z[1] + A[8] * z[2];
542: (void)PetscLogFlops(15.0);
543: }
544: /*
545: A: packed, row-major m x n array
546: x: length m array
547: y: length n arra
548: ldx: the stride in x and y
550: A[i,j] = A[i * n + j]
551: A^T[j,i] = A[i,j]
552: */
553: static inline void DMPlex_MultTransposeReal_Internal(const PetscReal A[], PetscInt m, PetscInt n, PetscInt ldx, const PetscScalar x[], PetscScalar y[])
554: {
555: PetscScalar z[3];
556: PetscInt i, j;
557: for (i = 0; i < m; ++i) z[i] = x[i * ldx];
558: for (j = 0; j < n; ++j) {
559: const PetscInt l = j * ldx;
560: y[l] = 0;
561: for (i = 0; i < m; ++i) y[l] += A[i * n + j] * z[i];
562: }
563: (void)PetscLogFlops(2 * m * n);
564: }
565: static inline void DMPlex_MultTranspose2DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
566: {
567: const PetscScalar z[2] = {x[0 * ldx], x[1 * ldx]};
568: y[0 * ldx] = A[0] * z[0] + A[2] * z[1];
569: y[1 * ldx] = A[1] * z[0] + A[3] * z[1];
570: (void)PetscLogFlops(6.0);
571: }
572: static inline void DMPlex_MultTranspose3DReal_Internal(const PetscReal A[], PetscInt ldx, const PetscScalar x[], PetscScalar y[])
573: {
574: const PetscScalar z[3] = {x[0 * ldx], x[1 * ldx], x[2 * ldx]};
575: y[0 * ldx] = A[0] * z[0] + A[3] * z[1] + A[6] * z[2];
576: y[1 * ldx] = A[1] * z[0] + A[4] * z[1] + A[7] * z[2];
577: y[2 * ldx] = A[2] * z[0] + A[5] * z[1] + A[8] * z[2];
578: (void)PetscLogFlops(15.0);
579: }
581: static inline void DMPlex_MatMult2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
582: {
583: #define PLEX_DIM__ 2
584: PetscScalar z[PLEX_DIM__];
585: for (PetscInt j = 0; j < n; ++j) {
586: for (int d = 0; d < PLEX_DIM__; ++d) z[d] = B[d * ldb + j];
587: DMPlex_Mult2D_Internal(A, 1, z, z);
588: for (int d = 0; d < PLEX_DIM__; ++d) C[d * ldb + j] = z[d];
589: }
590: (void)PetscLogFlops(8.0 * n);
591: #undef PLEX_DIM__
592: }
593: static inline void DMPlex_MatMult3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
594: {
595: #define PLEX_DIM__ 3
596: PetscScalar z[PLEX_DIM__];
597: for (PetscInt j = 0; j < n; ++j) {
598: for (int d = 0; d < PLEX_DIM__; ++d) z[d] = B[d * ldb + j];
599: DMPlex_Mult3D_Internal(A, 1, z, z);
600: for (int d = 0; d < PLEX_DIM__; ++d) C[d * ldb + j] = z[d];
601: }
602: (void)PetscLogFlops(8.0 * n);
603: #undef PLEX_DIM__
604: }
605: static inline void DMPlex_MatMultTranspose2D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
606: {
607: #define PLEX_DIM__ 2
608: PetscScalar z[PLEX_DIM__];
609: for (PetscInt j = 0; j < n; ++j) {
610: for (int d = 0; d < PLEX_DIM__; ++d) z[d] = B[d * ldb + j];
611: DMPlex_MultTranspose2D_Internal(A, 1, z, z);
612: for (int d = 0; d < PLEX_DIM__; ++d) C[d * ldb + j] = z[d];
613: }
614: (void)PetscLogFlops(8.0 * n);
615: #undef PLEX_DIM__
616: }
617: static inline void DMPlex_MatMultTranspose3D_Internal(const PetscScalar A[], PetscInt n, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
618: {
619: #define PLEX_DIM__ 3
620: PetscScalar z[PLEX_DIM__];
621: for (PetscInt j = 0; j < n; ++j) {
622: for (int d = 0; d < PLEX_DIM__; ++d) z[d] = B[d * ldb + j];
623: DMPlex_MultTranspose3D_Internal(A, 1, z, z);
624: for (int d = 0; d < PLEX_DIM__; ++d) C[d * ldb + j] = z[d];
625: }
626: (void)PetscLogFlops(8.0 * n);
627: #undef PLEX_DIM__
628: }
630: static inline void DMPlex_MatMultLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
631: {
632: for (PetscInt j = 0; j < m; ++j) DMPlex_MultTranspose2D_Internal(A, 1, &B[j * ldb], &C[j * ldb]);
633: (void)PetscLogFlops(8.0 * m);
634: }
635: static inline void DMPlex_MatMultLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
636: {
637: for (PetscInt j = 0; j < m; ++j) DMPlex_MultTranspose3D_Internal(A, 1, &B[j * ldb], &C[j * ldb]);
638: (void)PetscLogFlops(8.0 * m);
639: }
640: static inline void DMPlex_MatMultTransposeLeft2D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
641: {
642: for (PetscInt j = 0; j < m; ++j) DMPlex_Mult2D_Internal(A, 1, &B[j * ldb], &C[j * ldb]);
643: (void)PetscLogFlops(8.0 * m);
644: }
645: static inline void DMPlex_MatMultTransposeLeft3D_Internal(const PetscScalar A[], PetscInt m, PetscInt ldb, const PetscScalar B[], PetscScalar C[])
646: {
647: for (PetscInt j = 0; j < m; ++j) DMPlex_Mult3D_Internal(A, 1, &B[j * ldb], &C[j * ldb]);
648: (void)PetscLogFlops(8.0 * m);
649: }
651: static inline void DMPlex_PTAP2DReal_Internal(const PetscReal P[], const PetscScalar A[], PetscScalar C[])
652: {
653: PetscScalar out[4];
654: PetscInt i, j, k, l;
655: for (i = 0; i < 2; ++i) {
656: for (j = 0; j < 2; ++j) {
657: out[i * 2 + j] = 0.;
658: for (k = 0; k < 2; ++k) {
659: for (l = 0; l < 2; ++l) out[i * 2 + j] += P[k * 2 + i] * A[k * 2 + l] * P[l * 2 + j];
660: }
661: }
662: }
663: for (i = 0; i < 2 * 2; ++i) C[i] = out[i];
664: (void)PetscLogFlops(48.0);
665: }
666: static inline void DMPlex_PTAP3DReal_Internal(const PetscReal P[], const PetscScalar A[], PetscScalar C[])
667: {
668: PetscScalar out[9];
669: PetscInt i, j, k, l;
670: for (i = 0; i < 3; ++i) {
671: for (j = 0; j < 3; ++j) {
672: out[i * 3 + j] = 0.;
673: for (k = 0; k < 3; ++k) {
674: for (l = 0; l < 3; ++l) out[i * 3 + j] += P[k * 3 + i] * A[k * 3 + l] * P[l * 3 + j];
675: }
676: }
677: }
678: for (i = 0; i < 2 * 2; ++i) C[i] = out[i];
679: (void)PetscLogFlops(243.0);
680: }
681: /* TODO Fix for aliasing of A and C */
682: static inline void DMPlex_PTAPReal_Internal(const PetscReal P[], PetscInt m, PetscInt n, const PetscScalar A[], PetscScalar C[])
683: {
684: PetscInt i, j, k, l;
685: for (i = 0; i < n; ++i) {
686: for (j = 0; j < n; ++j) {
687: C[i * n + j] = 0.;
688: for (k = 0; k < m; ++k) {
689: for (l = 0; l < m; ++l) C[i * n + j] += P[k * n + i] * A[k * m + l] * P[l * n + j];
690: }
691: }
692: }
693: (void)PetscLogFlops(243.0);
694: }
696: static inline void DMPlex_Transpose2D_Internal(PetscScalar A[])
697: {
698: PetscScalar tmp;
699: tmp = A[1];
700: A[1] = A[2];
701: A[2] = tmp;
702: }
703: static inline void DMPlex_Transpose3D_Internal(PetscScalar A[])
704: {
705: PetscScalar tmp;
706: tmp = A[1];
707: A[1] = A[3];
708: A[3] = tmp;
709: tmp = A[2];
710: A[2] = A[6];
711: A[6] = tmp;
712: tmp = A[5];
713: A[5] = A[7];
714: A[7] = tmp;
715: }
717: static inline void DMPlex_Invert2D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
718: {
719: // Allow zero volume cells
720: const PetscReal invDet = detJ == 0 ? 1.0 : 1.0 / detJ;
722: invJ[0] = invDet * J[3];
723: invJ[1] = -invDet * J[1];
724: invJ[2] = -invDet * J[2];
725: invJ[3] = invDet * J[0];
726: (void)PetscLogFlops(5.0);
727: }
729: static inline void DMPlex_Invert3D_Internal(PetscReal invJ[], PetscReal J[], PetscReal detJ)
730: {
731: // Allow zero volume cells
732: const PetscReal invDet = detJ == 0 ? 1.0 : 1.0 / detJ;
734: invJ[0 * 3 + 0] = invDet * (J[1 * 3 + 1] * J[2 * 3 + 2] - J[1 * 3 + 2] * J[2 * 3 + 1]);
735: invJ[0 * 3 + 1] = invDet * (J[0 * 3 + 2] * J[2 * 3 + 1] - J[0 * 3 + 1] * J[2 * 3 + 2]);
736: invJ[0 * 3 + 2] = invDet * (J[0 * 3 + 1] * J[1 * 3 + 2] - J[0 * 3 + 2] * J[1 * 3 + 1]);
737: invJ[1 * 3 + 0] = invDet * (J[1 * 3 + 2] * J[2 * 3 + 0] - J[1 * 3 + 0] * J[2 * 3 + 2]);
738: invJ[1 * 3 + 1] = invDet * (J[0 * 3 + 0] * J[2 * 3 + 2] - J[0 * 3 + 2] * J[2 * 3 + 0]);
739: invJ[1 * 3 + 2] = invDet * (J[0 * 3 + 2] * J[1 * 3 + 0] - J[0 * 3 + 0] * J[1 * 3 + 2]);
740: invJ[2 * 3 + 0] = invDet * (J[1 * 3 + 0] * J[2 * 3 + 1] - J[1 * 3 + 1] * J[2 * 3 + 0]);
741: invJ[2 * 3 + 1] = invDet * (J[0 * 3 + 1] * J[2 * 3 + 0] - J[0 * 3 + 0] * J[2 * 3 + 1]);
742: invJ[2 * 3 + 2] = invDet * (J[0 * 3 + 0] * J[1 * 3 + 1] - J[0 * 3 + 1] * J[1 * 3 + 0]);
743: (void)PetscLogFlops(37.0);
744: }
746: static inline void DMPlex_Det2D_Internal(PetscReal *detJ, const PetscReal J[])
747: {
748: *detJ = J[0] * J[3] - J[1] * J[2];
749: (void)PetscLogFlops(3.0);
750: }
752: static inline void DMPlex_Det3D_Internal(PetscReal *detJ, const PetscReal J[])
753: {
754: *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]));
755: (void)PetscLogFlops(12.0);
756: }
758: static inline void DMPlex_Det2D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
759: {
760: *detJ = PetscRealPart(J[0]) * PetscRealPart(J[3]) - PetscRealPart(J[1]) * PetscRealPart(J[2]);
761: (void)PetscLogFlops(3.0);
762: }
764: static inline void DMPlex_Det3D_Scalar_Internal(PetscReal *detJ, const PetscScalar J[])
765: {
766: *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])));
767: (void)PetscLogFlops(12.0);
768: }
770: static inline void DMPlex_WaxpyD_Internal(PetscInt dim, PetscReal a, const PetscReal *x, const PetscReal *y, PetscReal *w)
771: {
772: PetscInt d;
773: for (d = 0; d < dim; ++d) w[d] = a * x[d] + y[d];
774: }
776: static inline PetscReal DMPlex_DotD_Internal(PetscInt dim, const PetscScalar *x, const PetscReal *y)
777: {
778: PetscReal sum = 0.0;
779: PetscInt d;
780: for (d = 0; d < dim; ++d) sum += PetscRealPart(x[d]) * y[d];
781: return sum;
782: }
784: static inline PetscReal DMPlex_DotRealD_Internal(PetscInt dim, const PetscReal *x, const PetscReal *y)
785: {
786: PetscReal sum = 0.0;
787: PetscInt d;
788: for (d = 0; d < dim; ++d) sum += x[d] * y[d];
789: return sum;
790: }
792: static inline PetscReal DMPlex_NormD_Internal(PetscInt dim, const PetscReal *x)
793: {
794: PetscReal sum = 0.0;
795: PetscInt d;
796: for (d = 0; d < dim; ++d) sum += x[d] * x[d];
797: return PetscSqrtReal(sum);
798: }
800: static inline PetscReal DMPlex_DistD_Internal(PetscInt dim, const PetscScalar *x, const PetscScalar *y)
801: {
802: PetscReal sum = 0.0;
803: PetscInt d;
804: for (d = 0; d < dim; ++d) sum += PetscRealPart(PetscConj(x[d] - y[d]) * (x[d] - y[d]));
805: return PetscSqrtReal(sum);
806: }
808: PETSC_INTERN PetscErrorCode DMPlexGetPointDualSpaceFEM(DM, PetscInt, PetscInt, PetscDualSpace *);
809: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPoint_Internal(PetscSection, PetscBool, PetscInt, PetscInt, PetscInt *, PetscBool, const PetscInt[], const PetscInt[], PetscInt[]);
810: PETSC_INTERN PetscErrorCode DMPlexGetIndicesPointFields_Internal(PetscSection, PetscBool, PetscInt, PetscInt, PetscInt[], PetscBool, const PetscInt ***, PetscInt, const PetscInt[], PetscInt[]);
811: PETSC_INTERN PetscErrorCode DMPlexGetTransitiveClosure_Internal(DM, PetscInt, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
812: PETSC_INTERN PetscErrorCode DMPlexVecGetOrientedClosure_Internal(DM, PetscSection, PetscBool, Vec, PetscInt, PetscInt, PetscInt *, PetscScalar *[]);
813: PETSC_INTERN PetscErrorCode DMPlexMatSetClosure_Internal(DM, PetscSection, PetscSection, PetscBool, Mat, PetscInt, const PetscScalar[], InsertMode);
815: PETSC_EXTERN PetscErrorCode DMPlexGetAllCells_Internal(DM, IS *);
816: PETSC_EXTERN PetscErrorCode DMPlexGetAllFaces_Internal(DM, IS *);
817: PETSC_EXTERN PetscErrorCode DMSNESGetFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
818: PETSC_EXTERN PetscErrorCode DMSNESRestoreFEGeom(DMField, IS, PetscQuadrature, PetscBool, PetscFEGeom **);
819: PETSC_EXTERN PetscErrorCode DMPlexComputeResidual_Patch_Internal(DM, PetscSection, IS, PetscReal, Vec, Vec, Vec, void *);
820: PETSC_EXTERN PetscErrorCode DMPlexComputeJacobian_Patch_Internal(DM, PetscSection, PetscSection, IS, PetscReal, PetscReal, Vec, Vec, Mat, Mat, void *);
821: PETSC_INTERN PetscErrorCode DMCreateSubDomainDM_Plex(DM, DMLabel, PetscInt, IS *, DM *);
822: PETSC_INTERN PetscErrorCode DMPlexBasisTransformPoint_Internal(DM, DM, Vec, PetscInt, PetscBool[], PetscBool, PetscScalar *);
823: PETSC_EXTERN PetscErrorCode DMPlexBasisTransformPointTensor_Internal(DM, DM, Vec, PetscInt, PetscBool, PetscInt, PetscScalar *);
824: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApplyReal_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscReal *, PetscReal *, void *);
825: PETSC_INTERN PetscErrorCode DMPlexBasisTransformApply_Internal(DM, const PetscReal[], PetscBool, PetscInt, const PetscScalar *, PetscScalar *, void *);
826: PETSC_INTERN PetscErrorCode DMCreateNeumannOverlap_Plex(DM, IS *, Mat *, PetscErrorCode (**)(Mat, PetscReal, Vec, Vec, PetscReal, IS, void *), void **);
827: PETSC_INTERN PetscErrorCode DMPlexMarkBoundaryFaces_Internal(DM, PetscInt, PetscInt, DMLabel, PetscBool);
828: PETSC_INTERN PetscErrorCode DMPlexDistributeOverlap_Internal(DM, PetscInt, MPI_Comm, const char *, PetscSF *, DM *);
830: PETSC_INTERN PetscErrorCode DMPeriodicCoordinateSetUp_Internal(DM);
832: /* Functions in the vtable */
833: PETSC_INTERN PetscErrorCode DMCreateInterpolation_Plex(DM dmCoarse, DM dmFine, Mat *interpolation, Vec *scaling);
834: PETSC_INTERN PetscErrorCode DMCreateInjection_Plex(DM dmCoarse, DM dmFine, Mat *mat);
835: PETSC_INTERN PetscErrorCode DMCreateMassMatrix_Plex(DM dmCoarse, DM dmFine, Mat *mat);
836: PETSC_INTERN PetscErrorCode DMCreateMassMatrixLumped_Plex(DM dmCoarse, Vec *lm);
837: PETSC_INTERN PetscErrorCode DMCreateLocalSection_Plex(DM dm);
838: PETSC_INTERN PetscErrorCode DMCreateDefaultConstraints_Plex(DM dm);
839: PETSC_INTERN PetscErrorCode DMCreateMatrix_Plex(DM dm, Mat *J);
840: PETSC_INTERN PetscErrorCode DMCreateCoordinateDM_Plex(DM dm, DM *cdm);
841: PETSC_INTERN PetscErrorCode DMCreateCoordinateField_Plex(DM dm, DMField *field);
842: PETSC_INTERN PetscErrorCode DMClone_Plex(DM dm, DM *newdm);
843: PETSC_INTERN PetscErrorCode DMSetUp_Plex(DM dm);
844: PETSC_INTERN PetscErrorCode DMDestroy_Plex(DM dm);
845: PETSC_INTERN PetscErrorCode DMView_Plex(DM dm, PetscViewer viewer);
846: PETSC_INTERN PetscErrorCode DMLoad_Plex(DM dm, PetscViewer viewer);
847: PETSC_INTERN PetscErrorCode DMCreateSubDM_Plex(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm);
848: PETSC_INTERN PetscErrorCode DMCreateSuperDM_Plex(DM dms[], PetscInt len, IS **is, DM *superdm);
849: PETSC_INTERN PetscErrorCode DMCreateDomainDecompositionScatters_Plex(DM, PetscInt, DM *, VecScatter **, VecScatter **, VecScatter **);
850: PETSC_INTERN PetscErrorCode DMCreateDomainDecomposition_Plex(DM, PetscInt *, char ***, IS **, IS **, DM **);
851: PETSC_INTERN PetscErrorCode DMCreateSectionPermutation_Plex(DM dm, IS *permutation, PetscBT *blockStarts);
853: // Coordinate mapping functions
854: PETSC_INTERN void coordMap_identity(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[]);
855: PETSC_INTERN void coordMap_shear(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[]);
856: PETSC_INTERN void coordMap_flare(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[]);
857: PETSC_INTERN void coordMap_annulus(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[]);
858: PETSC_INTERN void coordMap_shell(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[]);