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