Actual source code: petscdmplex.h
1: /*
2: DMPlex, for parallel unstructured distributed mesh problems.
3: */
4: #ifndef PETSCDMPLEX_H
5: #define PETSCDMPLEX_H
7: #include <petscsection.h>
8: #include <petscpartitioner.h>
9: #include <petscdm.h>
10: #include <petscdmplextypes.h>
11: #include <petscdt.h>
12: #include <petscfe.h>
13: #include <petscfv.h>
14: #include <petscsftypes.h>
15: #include <petscdmfield.h>
16: #include <petscviewer.h>
18: /* SUBMANSEC = DMPlex */
20: PETSC_EXTERN PetscErrorCode PetscPartitionerDMPlexPartition(PetscPartitioner, DM, PetscSection, PetscSection, IS *);
22: PETSC_EXTERN PetscErrorCode DMPlexBuildFromCellList(DM, PetscInt, PetscInt, PetscInt, const PetscInt[]);
23: PETSC_EXTERN PetscErrorCode DMPlexBuildFromCellListParallel(DM, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscSF *, PetscInt **);
24: PETSC_EXTERN PetscErrorCode DMPlexBuildCoordinatesFromCellList(DM, PetscInt, const PetscReal[]);
25: PETSC_EXTERN PetscErrorCode DMPlexBuildCoordinatesFromCellListParallel(DM, PetscInt, PetscSF, const PetscReal[]);
27: PETSC_EXTERN PetscErrorCode DMPlexCreate(MPI_Comm, DM *);
28: PETSC_EXTERN PetscErrorCode DMPlexCreateCohesiveSubmesh(DM, PetscBool, const char[], PetscInt, DM *);
29: PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellListPetsc(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const PetscInt[], PetscInt, const PetscReal[], DM *);
30: PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellListParallelPetsc(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const PetscInt[], PetscInt, const PetscReal[], PetscSF *, PetscInt **, DM *);
31: PETSC_EXTERN PetscErrorCode DMPlexCreateFromDAG(DM, PetscInt, const PetscInt[], const PetscInt[], const PetscInt[], const PetscInt[], const PetscScalar[]);
32: PETSC_EXTERN PetscErrorCode DMPlexCreateReferenceCell(MPI_Comm, DMPolytopeType, DM *);
33: PETSC_EXTERN PetscErrorCode DMPlexSetOptionsPrefix(DM, const char[]);
34: PETSC_EXTERN PetscErrorCode DMPlexGetChart(DM, PetscInt *, PetscInt *);
35: PETSC_EXTERN PetscErrorCode DMPlexSetChart(DM, PetscInt, PetscInt);
36: PETSC_EXTERN PetscErrorCode DMPlexGetConeSize(DM, PetscInt, PetscInt *);
37: PETSC_EXTERN PetscErrorCode DMPlexSetConeSize(DM, PetscInt, PetscInt);
38: PETSC_EXTERN PetscErrorCode DMPlexGetCone(DM, PetscInt, const PetscInt *[]);
39: PETSC_EXTERN PetscErrorCode DMPlexGetConeTuple(DM, IS, PetscSection *, IS *);
40: PETSC_EXTERN PetscErrorCode DMPlexGetConeRecursive(DM, IS, PetscInt *, IS *[], PetscSection *[]);
41: PETSC_EXTERN PetscErrorCode DMPlexRestoreConeRecursive(DM, IS, PetscInt *, IS *[], PetscSection *[]);
42: PETSC_EXTERN PetscErrorCode DMPlexGetConeRecursiveVertices(DM, IS, IS *);
43: PETSC_EXTERN PetscErrorCode DMPlexGetOrientedCone(DM, PetscInt, const PetscInt *[], const PetscInt *[]);
44: PETSC_EXTERN PetscErrorCode DMPlexRestoreOrientedCone(DM, PetscInt, const PetscInt *[], const PetscInt *[]);
45: PETSC_EXTERN PetscErrorCode DMPlexSetCone(DM, PetscInt, const PetscInt[]);
46: PETSC_EXTERN PetscErrorCode DMPlexInsertCone(DM, PetscInt, PetscInt, PetscInt);
47: PETSC_EXTERN PetscErrorCode DMPlexInsertConeOrientation(DM, PetscInt, PetscInt, PetscInt);
48: PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientation(DM, PetscInt, const PetscInt *[]);
49: PETSC_EXTERN PetscErrorCode DMPlexSetConeOrientation(DM, PetscInt, const PetscInt[]);
50: PETSC_EXTERN PetscErrorCode DMPlexGetSupportSize(DM, PetscInt, PetscInt *);
51: PETSC_EXTERN PetscErrorCode DMPlexSetSupportSize(DM, PetscInt, PetscInt);
52: PETSC_EXTERN PetscErrorCode DMPlexGetSupport(DM, PetscInt, const PetscInt *[]);
53: PETSC_EXTERN PetscErrorCode DMPlexSetSupport(DM, PetscInt, const PetscInt[]);
54: PETSC_EXTERN PetscErrorCode DMPlexInsertSupport(DM, PetscInt, PetscInt, PetscInt);
55: PETSC_EXTERN PetscErrorCode DMPlexGetConeSection(DM, PetscSection *);
56: PETSC_EXTERN PetscErrorCode DMPlexGetSupportSection(DM, PetscSection *);
57: PETSC_EXTERN PetscErrorCode DMPlexGetCones(DM, PetscInt *[]);
58: PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientations(DM, PetscInt *[]);
59: PETSC_EXTERN PetscErrorCode DMPlexGetMaxSizes(DM, PetscInt *, PetscInt *);
60: PETSC_EXTERN PetscErrorCode DMPlexSymmetrize(DM);
61: PETSC_EXTERN PetscErrorCode DMPlexStratify(DM);
62: PETSC_EXTERN PetscErrorCode DMPlexEqual(DM, DM, PetscBool *);
63: PETSC_EXTERN PetscErrorCode DMPlexOrientPoint(DM, PetscInt, PetscInt);
64: PETSC_EXTERN PetscErrorCode DMPlexOrient(DM);
65: PETSC_EXTERN PetscErrorCode DMPlexPreallocateOperator(DM, PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscInt[], Mat, PetscBool);
66: PETSC_EXTERN PetscErrorCode DMPlexGetPointLocal(DM, PetscInt, PetscInt *, PetscInt *);
67: PETSC_EXTERN PetscErrorCode DMPlexPointLocalRead(DM, PetscInt, const PetscScalar *, void *);
68: PETSC_EXTERN PetscErrorCode DMPlexPointLocalRef(DM, PetscInt, PetscScalar *, void *);
69: PETSC_EXTERN PetscErrorCode DMPlexGetPointLocalField(DM, PetscInt, PetscInt, PetscInt *, PetscInt *);
70: PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRef(DM, PetscInt, PetscInt, PetscScalar *, void *);
71: PETSC_EXTERN PetscErrorCode DMPlexPointLocalFieldRead(DM, PetscInt, PetscInt, const PetscScalar *, void *);
72: PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobal(DM, PetscInt, PetscInt *, PetscInt *);
73: PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRead(DM, PetscInt, const PetscScalar *, const void *);
74: PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRef(DM, PetscInt, PetscScalar *, void *);
75: PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobalField(DM, PetscInt, PetscInt, PetscInt *, PetscInt *);
76: PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRef(DM, PetscInt, PetscInt, PetscScalar *, void *);
77: PETSC_EXTERN PetscErrorCode DMPlexPointGlobalFieldRead(DM, PetscInt, PetscInt, const PetscScalar *, void *);
79: /* Topological interpolation */
80: PETSC_EXTERN const char *const DMPlexInterpolatedFlags[];
82: /*E
83: DMPlexInterpolatedFlag - Describes level of topological interpolatedness.
85: Local or collective property depending on whether it is returned by `DMPlexIsInterpolated()` or `DMPlexIsInterpolatedCollective()`.
87: $ DMPLEX_INTERPOLATED_INVALID - Uninitialized value (internal use only; never returned by `DMPlexIsInterpolated()` or `DMPlexIsInterpolatedCollective()`)
88: $ DMPLEX_INTERPOLATED_NONE - Mesh is not interpolated
89: $ DMPLEX_INTERPOLATED_PARTIAL - Mesh is partially interpolated. This can e.g. mean `DMPLEX` with cells, faces and vertices but no edges represented, or a mesh with mixed cones (see `DMPlexStratify()` for an example)
90: $ DMPLEX_INTERPOLATED_MIXED - Can be returned only by `DMPlexIsInterpolatedCollective()`, meaning that `DMPlexIsInterpolated()` returns different interpolatedness on different ranks
91: $ DMPLEX_INTERPOLATED_FULL - Mesh is fully interpolated
93: Level: intermediate
95: Note:
96: An interpolated `DMPLEX` means that edges (and faces for 3d meshes) are present in the `DMPLEX` data structures.
98: Developer Note:
99: Any additions/changes here MUST also be made in include/petsc/finclude/petscdmplex.h and src/dm/f90-mod/petscdmplex.h
101: .seealso: `DMPLEX`, `DMPlexIsInterpolated()`, `DMPlexIsInterpolatedCollective()`, `DMPlexInterpolate()`, `DMPlexUninterpolate()`
102: E*/
103: typedef enum {
104: DMPLEX_INTERPOLATED_INVALID = -1,
105: DMPLEX_INTERPOLATED_NONE = 0,
106: DMPLEX_INTERPOLATED_PARTIAL,
107: DMPLEX_INTERPOLATED_MIXED,
108: DMPLEX_INTERPOLATED_FULL
109: } DMPlexInterpolatedFlag;
111: PETSC_EXTERN PetscErrorCode DMPlexInterpolate(DM, DM *);
112: PETSC_EXTERN PetscErrorCode DMPlexUninterpolate(DM, DM *);
113: PETSC_EXTERN PetscErrorCode DMPlexInterpolatePointSF(DM, PetscSF);
114: PETSC_EXTERN PetscErrorCode DMPlexIsInterpolated(DM, DMPlexInterpolatedFlag *);
115: PETSC_EXTERN PetscErrorCode DMPlexIsInterpolatedCollective(DM, DMPlexInterpolatedFlag *);
117: PETSC_EXTERN PetscErrorCode DMPlexFilter(DM, DMLabel, PetscInt, DM *);
118: PETSC_EXTERN PetscErrorCode DMPlexGetCellNumbering(DM, IS *);
119: PETSC_EXTERN PetscErrorCode DMPlexGetVertexNumbering(DM, IS *);
120: PETSC_EXTERN PetscErrorCode DMPlexCreatePointNumbering(DM, IS *);
121: PETSC_EXTERN PetscErrorCode DMPlexCreateRankField(DM, Vec *);
122: PETSC_EXTERN PetscErrorCode DMPlexCreateLabelField(DM, DMLabel, Vec *);
124: PETSC_EXTERN PetscErrorCode DMPlexGetDepth(DM, PetscInt *);
125: PETSC_EXTERN PetscErrorCode DMPlexGetDepthLabel(DM, DMLabel *);
126: PETSC_EXTERN PetscErrorCode DMPlexGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *);
127: PETSC_EXTERN PetscErrorCode DMPlexGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
128: PETSC_EXTERN PetscErrorCode DMPlexGetPointDepth(DM, PetscInt, PetscInt *);
129: PETSC_EXTERN PetscErrorCode DMPlexGetPointHeight(DM, PetscInt, PetscInt *);
130: PETSC_EXTERN PetscErrorCode DMPlexGetCellTypeLabel(DM, DMLabel *);
131: PETSC_EXTERN PetscErrorCode DMPlexGetCellType(DM, PetscInt, DMPolytopeType *);
132: PETSC_EXTERN PetscErrorCode DMPlexSetCellType(DM, PetscInt, DMPolytopeType);
133: PETSC_EXTERN PetscErrorCode DMPlexComputeCellTypes(DM);
134: PETSC_EXTERN PetscErrorCode DMPlexInvertCell(DMPolytopeType, PetscInt[]);
135: PETSC_EXTERN PetscErrorCode DMPlexReorderCell(DM, PetscInt, PetscInt[]);
137: /* Topological Operations */
138: PETSC_EXTERN PetscErrorCode DMPlexGetMeet(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **);
139: PETSC_EXTERN PetscErrorCode DMPlexGetFullMeet(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **);
140: PETSC_EXTERN PetscErrorCode DMPlexRestoreMeet(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **);
141: PETSC_EXTERN PetscErrorCode DMPlexGetJoin(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **);
142: PETSC_EXTERN PetscErrorCode DMPlexGetFullJoin(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **);
143: PETSC_EXTERN PetscErrorCode DMPlexRestoreJoin(DM, PetscInt, const PetscInt[], PetscInt *, const PetscInt **);
144: PETSC_EXTERN PetscErrorCode DMPlexGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
145: PETSC_EXTERN PetscErrorCode DMPlexRestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
146: PETSC_EXTERN PetscErrorCode DMPlexGetCompressedClosure(DM, PetscSection, PetscInt, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);
147: PETSC_EXTERN PetscErrorCode DMPlexRestoreCompressedClosure(DM, PetscSection, PetscInt, PetscInt *, PetscInt **, PetscSection *, IS *, const PetscInt **);
149: /*E
150: DMPlexTPSType - Type of triply-periodic surface for a `DMPLEX`
152: $ DMPLEX_TPS_SCHWARZ_P - Schwarz Primitive surface, defined by the equation cos(x) + cos(y) + cos(z) = 0.
153: $ DMPLEX_TPS_GYROID - Gyroid surface, defined by the equation sin(x)cos(y) + sin(y)cos(z) + sin(z)cos(x) = 0
155: Level: intermediate
157: Developer Note:
158: Any additions/changes here MUST also be made in include/petsc/finclude/petscdmplex.h and src/dm/f90-mod/petscdmplex.h
160: .seealso: `DMPlexCreateTPSMesh()`
161: E*/
162: typedef enum {
163: DMPLEX_TPS_SCHWARZ_P,
164: DMPLEX_TPS_GYROID
165: } DMPlexTPSType;
166: PETSC_EXTERN const char *const DMPlexTPSTypes[];
168: /* Mesh Generation */
169: PETSC_EXTERN PetscErrorCode DMPlexGenerate(DM, const char[], PetscBool, DM *);
170: PETSC_EXTERN PetscErrorCode DMPlexCopyCoordinates(DM, DM);
171: PETSC_EXTERN PetscErrorCode DMPlexCreateCoordinateSpace(DM, PetscInt, 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[]));
172: PETSC_EXTERN PetscErrorCode DMPlexCreateDoublet(MPI_Comm, PetscInt, PetscBool, PetscBool, PetscReal, DM *);
173: PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, DM *);
174: PETSC_EXTERN PetscErrorCode DMPlexCreateBoxSurfaceMesh(MPI_Comm, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscBool, DM *);
175: PETSC_EXTERN PetscErrorCode DMPlexCreateSphereMesh(MPI_Comm, PetscInt, PetscBool, PetscReal, DM *);
176: PETSC_EXTERN PetscErrorCode DMPlexCreateBallMesh(MPI_Comm, PetscInt, PetscReal, DM *);
177: PETSC_EXTERN PetscErrorCode DMPlexCreateHexCylinderMesh(MPI_Comm, DMBoundaryType, DM *);
178: PETSC_EXTERN PetscErrorCode DMPlexCreateTPSMesh(MPI_Comm, DMPlexTPSType, const PetscInt[], const DMBoundaryType[], PetscBool, PetscInt, PetscInt, PetscReal, DM *);
179: PETSC_EXTERN PetscErrorCode DMPlexCreateWedgeCylinderMesh(MPI_Comm, PetscInt, PetscBool, DM *);
180: PETSC_EXTERN PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm, const PetscInt[], const PetscReal[], const PetscReal[], const DMBoundaryType[], PetscBool, PetscBool, DM *);
181: PETSC_EXTERN PetscErrorCode DMPlexExtrude(DM, PetscInt, PetscReal, PetscBool, PetscBool, const PetscReal[], const PetscReal[], DM *);
182: PETSC_EXTERN PetscErrorCode DMPlexCreateConeSection(DM, PetscSection *);
183: PETSC_EXTERN PetscErrorCode DMPlexInflateToGeomModel(DM);
185: PETSC_EXTERN PetscErrorCode DMPlexSetIsoperiodicFaceSF(DM, PetscSF);
186: PETSC_EXTERN PetscErrorCode DMPlexGetIsoperiodicFaceSF(DM, PetscSF *);
187: PETSC_EXTERN PetscErrorCode DMPlexSetIsoperiodicFaceTransform(DM, const PetscScalar[]);
189: PETSC_EXTERN PetscErrorCode DMPlexCheck(DM);
190: PETSC_EXTERN PetscErrorCode DMPlexCheckSymmetry(DM);
191: PETSC_EXTERN PetscErrorCode DMPlexCheckSkeleton(DM, PetscInt);
192: PETSC_EXTERN PetscErrorCode DMPlexCheckFaces(DM, PetscInt);
193: PETSC_EXTERN PetscErrorCode DMPlexCheckGeometry(DM);
194: PETSC_EXTERN PetscErrorCode DMPlexCheckPointSF(DM, PetscSF, PetscBool);
195: PETSC_EXTERN PetscErrorCode DMPlexCheckInterfaceCones(DM);
196: PETSC_EXTERN PetscErrorCode DMPlexCheckCellShape(DM, PetscBool, PetscReal);
197: PETSC_EXTERN PetscErrorCode DMPlexComputeOrthogonalQuality(DM, PetscFV, PetscReal, Vec *, DMLabel *);
199: PETSC_EXTERN PetscErrorCode DMPlexTriangleSetOptions(DM, const char *);
200: PETSC_EXTERN PetscErrorCode DMPlexTetgenSetOptions(DM, const char *);
202: PETSC_EXTERN PetscErrorCode DMPlexCreateFromFile(MPI_Comm, const char[], const char[], PetscBool, DM *);
203: PETSC_EXTERN PetscErrorCode DMPlexCreateExodus(MPI_Comm, PetscInt, PetscBool, DM *);
204: PETSC_EXTERN PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm, const char[], PetscBool, DM *);
205: PETSC_EXTERN PetscErrorCode DMPlexCreateCGNS(MPI_Comm, PetscInt, PetscBool, DM *);
206: PETSC_EXTERN PetscErrorCode DMPlexCreateCGNSFromFile(MPI_Comm, const char[], PetscBool, DM *);
207: PETSC_EXTERN PetscErrorCode DMPlexCreateGmsh(MPI_Comm, PetscViewer, PetscBool, DM *);
208: PETSC_EXTERN PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm, const char[], PetscBool, DM *);
209: PETSC_EXTERN PetscErrorCode DMPlexCreateFluent(MPI_Comm, PetscViewer, PetscBool, DM *);
210: PETSC_EXTERN PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm, const char[], PetscBool, DM *);
211: PETSC_EXTERN PetscErrorCode DMPlexCreateMedFromFile(MPI_Comm, const char[], PetscBool, DM *);
212: PETSC_EXTERN PetscErrorCode DMPlexCreatePLYFromFile(MPI_Comm, const char[], PetscBool, DM *);
213: PETSC_EXTERN PetscErrorCode DMPlexCreateEGADSFromFile(MPI_Comm, const char[], DM *);
214: PETSC_EXTERN PetscErrorCode DMPlexCreateEGADSLiteFromFile(MPI_Comm, const char[], DM *);
216: PETSC_EXTERN PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo);
217: PETSC_EXTERN PetscErrorCode PetscViewerExodusIIGetId(PetscViewer, int *);
218: PETSC_EXTERN PetscErrorCode PetscViewerExodusIISetOrder(PetscViewer, PetscInt);
219: PETSC_EXTERN PetscErrorCode PetscViewerExodusIIGetOrder(PetscViewer, PetscInt *);
221: /* Mesh Partitioning and Distribution */
222: #define DMPLEX_OVERLAP_MANUAL -1
224: PETSC_EXTERN PetscErrorCode DMPlexCreateNeighborCSR(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **);
225: PETSC_EXTERN PetscErrorCode DMPlexGetPartitioner(DM, PetscPartitioner *);
226: PETSC_EXTERN PetscErrorCode DMPlexSetPartitioner(DM, PetscPartitioner);
227: PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionerGraph(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **, IS *);
228: PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelInvert(DM, DMLabel, PetscSF, DMLabel);
229: PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelClosure(DM, DMLabel);
230: PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelAdjacency(DM, DMLabel);
231: PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelPropagate(DM, DMLabel);
232: PETSC_EXTERN PetscErrorCode DMPlexPartitionLabelCreateSF(DM, DMLabel, PetscSF *);
233: PETSC_EXTERN PetscErrorCode DMPlexSetPartitionBalance(DM, PetscBool);
234: PETSC_EXTERN PetscErrorCode DMPlexGetPartitionBalance(DM, PetscBool *);
235: PETSC_EXTERN PetscErrorCode DMPlexIsDistributed(DM, PetscBool *);
236: PETSC_EXTERN PetscErrorCode DMPlexDistribute(DM, PetscInt, PetscSF *, DM *);
237: PETSC_EXTERN PetscErrorCode DMPlexDistributeOverlap(DM, PetscInt, PetscSF *, DM *);
238: PETSC_EXTERN PetscErrorCode DMPlexGetOverlap(DM, PetscInt *);
239: PETSC_EXTERN PetscErrorCode DMPlexSetOverlap(DM, DM, PetscInt);
240: PETSC_EXTERN PetscErrorCode DMPlexDistributeGetDefault(DM, PetscBool *);
241: PETSC_EXTERN PetscErrorCode DMPlexDistributeSetDefault(DM, PetscBool);
242: PETSC_EXTERN PetscErrorCode DMPlexDistributeField(DM, PetscSF, PetscSection, Vec, PetscSection, Vec);
243: PETSC_EXTERN PetscErrorCode DMPlexDistributeFieldIS(DM, PetscSF, PetscSection, IS, PetscSection, IS *);
244: PETSC_EXTERN PetscErrorCode DMPlexDistributeData(DM, PetscSF, PetscSection, MPI_Datatype, void *, PetscSection, void **) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(5, 4) PETSC_ATTRIBUTE_MPI_POINTER_WITH_TYPE(7, 4);
245: PETSC_EXTERN PetscErrorCode DMPlexRebalanceSharedPoints(DM, PetscInt, PetscBool, PetscBool, PetscBool *);
246: PETSC_EXTERN PetscErrorCode DMPlexMigrate(DM, PetscSF, DM);
247: PETSC_EXTERN PetscErrorCode DMPlexGetGatherDM(DM, PetscSF *, DM *);
248: PETSC_EXTERN PetscErrorCode DMPlexGetRedundantDM(DM, PetscSF *, DM *);
249: PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUser(DM, PetscErrorCode (*)(DM, PetscInt, PetscInt *, PetscInt[], void *), void *);
250: PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUser(DM, PetscErrorCode (**)(DM, PetscInt, PetscInt *, PetscInt[], void *), void **);
251: PETSC_EXTERN PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM, PetscBool);
252: PETSC_EXTERN PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM, PetscBool *);
253: PETSC_EXTERN PetscErrorCode DMPlexGetAdjacency(DM, PetscInt, PetscInt *, PetscInt *[]);
254: PETSC_EXTERN PetscErrorCode DMPlexSetMigrationSF(DM, PetscSF);
255: PETSC_EXTERN PetscErrorCode DMPlexGetMigrationSF(DM, PetscSF *);
256: PETSC_EXTERN PetscErrorCode DMPlexDistributionSetName(DM, const char[]);
257: PETSC_EXTERN PetscErrorCode DMPlexDistributionGetName(DM, const char *[]);
259: /*E
260: DMPlexReorderDefaultFlag - Flag indicating whether the `DMPLEX` should be reordered by default
262: $ DMPLEX_REORDER_DEFAULT_NOTSET - Flag not set.
263: $ DMPLEX_REORDER_DEFAULT_FALSE - Do not reorder by default.
264: $ DMPLEX_REORDER_DEFAULT_TRUE - Reorder by default.
266: Level: intermediate
268: Developer Note:
269: Any additions/changes here MUST also be made in include/petsc/finclude/petscdmplex.h and src/dm/f90-mod/petscdmplex.h
271: Could be replaced with `PETSC_BOOL3`
273: .seealso: `DMPlexReorderSetDefault()`, `DMPlexReorderGetDefault()`, `DMPlexGetOrdering()`, `DMPlexPermute()`
274: E*/
275: typedef enum {
276: DMPLEX_REORDER_DEFAULT_NOTSET = -1,
277: DMPLEX_REORDER_DEFAULT_FALSE = 0,
278: DMPLEX_REORDER_DEFAULT_TRUE
279: } DMPlexReorderDefaultFlag;
280: PETSC_EXTERN PetscErrorCode DMPlexGetOrdering(DM, MatOrderingType, DMLabel, IS *);
281: PETSC_EXTERN PetscErrorCode DMPlexGetOrdering1D(DM, IS *);
282: PETSC_EXTERN PetscErrorCode DMPlexPermute(DM, IS, DM *);
283: PETSC_EXTERN PetscErrorCode DMPlexReorderGetDefault(DM, DMPlexReorderDefaultFlag *);
284: PETSC_EXTERN PetscErrorCode DMPlexReorderSetDefault(DM, DMPlexReorderDefaultFlag);
286: PETSC_EXTERN PetscErrorCode DMPlexCreateProcessSF(DM, PetscSF, IS *, PetscSF *);
287: PETSC_EXTERN PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM, PetscSF, PetscSection, IS, PetscSection, IS, IS *, PetscSF *);
288: PETSC_EXTERN PetscErrorCode DMPlexDistributeOwnership(DM, PetscSection, IS *, PetscSection, IS *);
289: PETSC_EXTERN PetscErrorCode DMPlexCreatePointSF(DM, PetscSF, PetscBool, PetscSF *);
290: PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapLabel(DM, PetscInt, PetscSection, IS, PetscSection, IS, DMLabel *);
291: PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapLabelFromLabels(DM, PetscInt, const DMLabel[], const PetscInt[], PetscInt, const DMLabel[], const PetscInt[], PetscSection, IS, PetscSection, IS, DMLabel *);
292: PETSC_EXTERN PetscErrorCode DMPlexCreateOverlapMigrationSF(DM, PetscSF, PetscSF *);
293: PETSC_EXTERN PetscErrorCode DMPlexStratifyMigrationSF(DM, PetscSF, PetscSF *);
295: /* Submesh Support */
296: PETSC_EXTERN PetscErrorCode DMPlexCreateSubmesh(DM, DMLabel, PetscInt, PetscBool, DM *);
297: PETSC_EXTERN PetscErrorCode DMPlexCreateHybridMesh(DM, DMLabel, DMLabel, PetscInt, DMLabel *, DMLabel *, DM *, DM *);
298: PETSC_EXTERN PetscErrorCode DMPlexGetSubpointMap(DM, DMLabel *);
299: PETSC_EXTERN PetscErrorCode DMPlexSetSubpointMap(DM, DMLabel);
300: PETSC_EXTERN PetscErrorCode DMPlexGetSubpointIS(DM, IS *);
302: PETSC_EXTERN PetscErrorCode DMGetEnclosureRelation(DM, DM, DMEnclosureType *);
303: PETSC_EXTERN PetscErrorCode DMGetEnclosurePoint(DM, DM, DMEnclosureType, PetscInt, PetscInt *);
305: PETSC_EXTERN PetscErrorCode DMPlexLabelComplete(DM, DMLabel);
306: PETSC_EXTERN PetscErrorCode DMPlexLabelCohesiveComplete(DM, DMLabel, DMLabel, PetscInt, PetscBool, DM);
307: PETSC_EXTERN PetscErrorCode DMPlexLabelAddCells(DM, DMLabel);
308: PETSC_EXTERN PetscErrorCode DMPlexLabelAddFaceCells(DM, DMLabel);
309: PETSC_EXTERN PetscErrorCode DMPlexLabelClearCells(DM, DMLabel);
311: PETSC_EXTERN PetscErrorCode DMPlexGetRefinementLimit(DM, PetscReal *);
312: PETSC_EXTERN PetscErrorCode DMPlexSetRefinementLimit(DM, PetscReal);
313: PETSC_EXTERN PetscErrorCode DMPlexGetRefinementUniform(DM, PetscBool *);
314: PETSC_EXTERN PetscErrorCode DMPlexSetRefinementUniform(DM, PetscBool);
315: PETSC_EXTERN PetscErrorCode DMPlexGetRefinementFunction(DM, PetscErrorCode (**)(const PetscReal[], PetscReal *));
316: PETSC_EXTERN PetscErrorCode DMPlexSetRefinementFunction(DM, PetscErrorCode (*)(const PetscReal[], PetscReal *));
317: PETSC_EXTERN PetscErrorCode DMPlexCreateCoarsePointIS(DM, IS *);
318: PETSC_EXTERN PetscErrorCode DMPlexGetRegularRefinement(DM, PetscBool *);
319: PETSC_EXTERN PetscErrorCode DMPlexSetRegularRefinement(DM, PetscBool);
321: /* Support for cell-vertex meshes */
322: PETSC_EXTERN PetscErrorCode DMPlexGetNumFaceVertices(DM, PetscInt, PetscInt, PetscInt *);
323: PETSC_EXTERN PetscErrorCode DMPlexGetOrientedFace(DM, PetscInt, PetscInt, const PetscInt[], PetscInt, PetscInt[], PetscInt[], PetscInt[], PetscBool *);
325: /* Geometry Support */
326: PETSC_EXTERN PetscErrorCode DMPlexGetMinRadius(DM, PetscReal *);
327: PETSC_EXTERN PetscErrorCode DMPlexSetMinRadius(DM, PetscReal);
328: PETSC_EXTERN PetscErrorCode DMPlexComputeProjection2Dto1D(PetscScalar[], PetscReal[]);
329: PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto1D(PetscScalar[], PetscReal[]);
330: PETSC_EXTERN PetscErrorCode DMPlexComputeProjection3Dto2D(PetscInt, PetscScalar[], PetscReal[]);
332: /* Point Location */
333: typedef struct _PetscGridHash *PetscGridHash;
334: PETSC_EXTERN PetscErrorCode PetscGridHashCreate(MPI_Comm, PetscInt, const PetscScalar[], PetscGridHash *);
335: PETSC_EXTERN PetscErrorCode PetscGridHashEnlarge(PetscGridHash, const PetscScalar[]);
336: PETSC_EXTERN PetscErrorCode PetscGridHashSetGrid(PetscGridHash, const PetscInt[], const PetscReal[]);
337: PETSC_EXTERN PetscErrorCode PetscGridHashGetEnclosingBox(PetscGridHash, PetscInt, const PetscScalar[], PetscInt[], PetscInt[]);
338: PETSC_EXTERN PetscErrorCode PetscGridHashDestroy(PetscGridHash *);
339: PETSC_EXTERN PetscErrorCode DMPlexFindVertices(DM, Vec, PetscReal, IS *);
341: /* FVM Support */
342: PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFVM(DM, PetscInt, PetscReal *, PetscReal[], PetscReal[]);
343: PETSC_EXTERN PetscErrorCode DMPlexComputeGeometryFVM(DM, Vec *, Vec *);
344: PETSC_EXTERN PetscErrorCode DMPlexComputeGradientFVM(DM, PetscFV, Vec, Vec, DM *);
345: PETSC_EXTERN PetscErrorCode DMPlexGetDataFVM(DM, PetscFV, Vec *, Vec *, DM *);
347: /* FEM Support */
348: PETSC_EXTERN PetscErrorCode DMPlexGetGeometryFVM(DM, Vec *, Vec *, PetscReal *);
349: PETSC_EXTERN PetscErrorCode DMPlexGetGradientDM(DM, PetscFV, DM *);
350: PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValues(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
351: PETSC_EXTERN PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues(DM, PetscBool, Vec, PetscReal, Vec, Vec, Vec);
352: PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssential(DM, PetscReal, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt[], PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, Vec);
353: PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssentialField(DM, PetscReal, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt[], 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[]), void *, Vec);
354: PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesEssentialBdField(DM, PetscReal, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt[], 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[]), void *, Vec);
355: PETSC_EXTERN PetscErrorCode DMPlexInsertBoundaryValuesRiemann(DM, PetscReal, Vec, Vec, Vec, PetscInt, PetscInt, const PetscInt[], DMLabel, PetscInt, const PetscInt[], PetscErrorCode (*)(PetscReal, const PetscReal *, const PetscReal *, const PetscScalar *, PetscScalar *, void *), void *, Vec);
356: PETSC_EXTERN PetscErrorCode DMPlexMarkBoundaryFaces(DM, PetscInt, DMLabel);
358: PETSC_EXTERN PetscErrorCode DMPlexCreateSection(DM, DMLabel[], const PetscInt[], const PetscInt[], PetscInt, const PetscInt[], const IS[], const IS[], IS, PetscSection *);
359: PETSC_EXTERN PetscErrorCode DMPlexGetSubdomainSection(DM, PetscSection *);
361: PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryAffineFEM(DM, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
362: PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometryFEM(DM, PetscInt, PetscQuadrature, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
363: PETSC_EXTERN PetscErrorCode DMPlexGetCellCoordinates(DM, PetscInt, PetscBool *, PetscInt *, const PetscScalar *[], PetscScalar *[]);
364: PETSC_EXTERN PetscErrorCode DMPlexRestoreCellCoordinates(DM, PetscInt, PetscBool *, PetscInt *, const PetscScalar *[], PetscScalar *[]);
365: PETSC_EXTERN PetscErrorCode DMPlexCoordinatesToReference(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[]);
366: PETSC_EXTERN PetscErrorCode DMPlexReferenceToCoordinates(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[]);
367: PETSC_EXTERN PetscErrorCode DMPlexShearGeometry(DM, DMDirection, PetscReal[]);
368: PETSC_EXTERN PetscErrorCode DMPlexRemapGeometry(DM, PetscReal, 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[]));
370: PETSC_EXTERN PetscErrorCode DMPlexVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]);
371: PETSC_EXTERN PetscErrorCode DMPlexVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]);
372: PETSC_EXTERN PetscErrorCode DMPlexVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar[], InsertMode);
373: PETSC_EXTERN PetscErrorCode DMPlexMatSetClosure(DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
374: PETSC_EXTERN PetscErrorCode DMPlexMatSetClosureGeneral(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
375: PETSC_EXTERN PetscErrorCode DMPlexGetClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscBool, PetscInt *, PetscInt *[], PetscInt *, PetscScalar *[]);
376: PETSC_EXTERN PetscErrorCode DMPlexRestoreClosureIndices(DM, PetscSection, PetscSection, PetscInt, PetscBool, PetscInt *, PetscInt *[], PetscInt *, PetscScalar *[]);
377: PETSC_EXTERN PetscErrorCode DMPlexMatSetClosureRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
378: PETSC_EXTERN PetscErrorCode DMPlexMatGetClosureIndicesRefined(DM, PetscSection, PetscSection, DM, PetscSection, PetscSection, PetscInt, PetscInt[], PetscInt[]);
379: PETSC_EXTERN PetscErrorCode DMPlexCreateClosureIndex(DM, PetscSection);
380: PETSC_EXTERN PetscErrorCode DMPlexSetClosurePermutationTensor(DM, PetscInt, PetscSection);
382: PETSC_EXTERN PetscErrorCode DMPlexConstructGhostCells(DM, const char[], PetscInt *, DM *);
383: PETSC_EXTERN PetscErrorCode DMPlexConstructCohesiveCells(DM, DMLabel, DMLabel, DM *);
384: PETSC_EXTERN PetscErrorCode DMPlexReorderCohesiveSupports(DM);
386: PETSC_EXTERN PetscErrorCode DMPlexGetVTKCellHeight(DM, PetscInt *);
387: PETSC_EXTERN PetscErrorCode DMPlexSetVTKCellHeight(DM, PetscInt);
388: PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll(PetscObject, PetscViewer);
389: PETSC_EXTERN PetscErrorCode DMPlexGetGhostCellStratum(DM, PetscInt *, PetscInt *);
390: PETSC_EXTERN PetscErrorCode DMPlexGetSimplexOrBoxCells(DM, PetscInt, PetscInt *, PetscInt *);
391: PETSC_EXTERN PetscErrorCode DMPlexIsSimplex(DM, PetscBool *);
393: PETSC_EXTERN PetscErrorCode DMPlexGetCellFields(DM, IS, Vec, Vec, Vec, PetscScalar **, PetscScalar **, PetscScalar **);
394: PETSC_EXTERN PetscErrorCode DMPlexRestoreCellFields(DM, IS, Vec, Vec, Vec, PetscScalar **, PetscScalar **, PetscScalar **);
395: PETSC_EXTERN PetscErrorCode DMPlexGetFaceFields(DM, PetscInt, PetscInt, Vec, Vec, Vec, Vec, Vec, PetscInt *, PetscScalar **, PetscScalar **);
396: PETSC_EXTERN PetscErrorCode DMPlexRestoreFaceFields(DM, PetscInt, PetscInt, Vec, Vec, Vec, Vec, Vec, PetscInt *, PetscScalar **, PetscScalar **);
397: PETSC_EXTERN PetscErrorCode DMPlexGetFaceGeometry(DM, PetscInt, PetscInt, Vec, Vec, PetscInt *, PetscFVFaceGeom **, PetscReal **);
398: PETSC_EXTERN PetscErrorCode DMPlexRestoreFaceGeometry(DM, PetscInt, PetscInt, Vec, Vec, PetscInt *, PetscFVFaceGeom **, PetscReal **);
400: PETSC_EXTERN PetscErrorCode DMPlexGetScale(DM, PetscUnit, PetscReal *);
401: PETSC_EXTERN PetscErrorCode DMPlexSetScale(DM, PetscUnit, PetscReal);
403: typedef struct {
404: DM dm;
405: Vec u; /* The base vector for the Jacobian action J(u) x */
406: Mat J; /* Preconditioner for testing */
407: void *user;
408: } JacActionCtx;
410: PETSC_EXTERN PetscErrorCode DMPlexSetMaxProjectionHeight(DM, PetscInt);
411: PETSC_EXTERN PetscErrorCode DMPlexGetMaxProjectionHeight(DM, PetscInt *);
412: PETSC_EXTERN PetscErrorCode DMPlexGetActivePoint(DM, PetscInt *);
413: PETSC_EXTERN PetscErrorCode DMPlexSetActivePoint(DM, PetscInt);
414: PETSC_EXTERN PetscErrorCode DMPlexProjectFieldLocal(DM, 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[], PetscScalar[]), InsertMode, Vec);
415: PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffLocal(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
416: PETSC_EXTERN PetscErrorCode DMPlexComputeL2FieldDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal[]);
417: PETSC_EXTERN PetscErrorCode DMPlexComputeL2DiffVec(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, Vec);
418: PETSC_EXTERN PetscErrorCode DMPlexComputeCellwiseIntegralFEM(DM, Vec, Vec, void *);
419: PETSC_EXTERN PetscErrorCode DMPlexComputeIntegralFEM(DM, Vec, PetscScalar *, void *);
420: PETSC_EXTERN PetscErrorCode DMPlexComputeBdIntegral(DM, Vec, DMLabel, PetscInt, const PetscInt[], 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[]), PetscScalar *, void *);
421: PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorNested(DM, DM, PetscBool, Mat, void *);
422: PETSC_EXTERN PetscErrorCode DMPlexComputeInterpolatorGeneral(DM, DM, Mat, void *);
423: PETSC_EXTERN PetscErrorCode DMPlexComputeClementInterpolant(DM, Vec, Vec);
424: PETSC_EXTERN PetscErrorCode DMPlexComputeGradientClementInterpolant(DM, Vec, Vec);
425: PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorFEM(DM, DM, VecScatter *, void *);
426: PETSC_EXTERN PetscErrorCode DMPlexComputeMassMatrixNested(DM, DM, Mat, void *);
427: PETSC_EXTERN PetscErrorCode DMPlexComputeMassMatrixGeneral(DM, DM, Mat, void *);
429: PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBody(DM, PetscInt, MatNullSpace *);
430: PETSC_EXTERN PetscErrorCode DMPlexCreateRigidBodies(DM, PetscInt, DMLabel, const PetscInt[], const PetscInt[], MatNullSpace *);
432: PETSC_EXTERN PetscErrorCode DMPlexSetSNESLocalFEM(DM, void *, void *, void *);
433: PETSC_EXTERN PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM, Vec, void *);
434: PETSC_EXTERN PetscErrorCode DMPlexSNESComputeResidualFEM(DM, Vec, Vec, void *);
435: PETSC_EXTERN PetscErrorCode DMPlexSNESComputeJacobianFEM(DM, Vec, Mat, Mat, void *);
436: PETSC_EXTERN PetscErrorCode DMPlexComputeBdResidualSingle(DM, PetscReal, PetscWeakForm, PetscFormKey, Vec, Vec, Vec);
437: PETSC_EXTERN PetscErrorCode DMPlexComputeBdJacobianSingle(DM, PetscReal, PetscWeakForm, DMLabel, PetscInt, const PetscInt[], PetscInt, Vec, Vec, PetscReal, Mat, Mat);
439: PETSC_EXTERN PetscErrorCode DMPlexTSComputeBoundary(DM, PetscReal, Vec, Vec, void *);
440: PETSC_EXTERN PetscErrorCode DMPlexTSComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *);
441: PETSC_EXTERN PetscErrorCode DMPlexTSComputeIFunctionFEM(DM, PetscReal, Vec, Vec, Vec, void *);
442: PETSC_EXTERN PetscErrorCode DMPlexTSComputeIJacobianFEM(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
443: PETSC_EXTERN PetscErrorCode DMPlexTSComputeRHSFunctionFEM(DM, PetscReal, Vec, Vec, void *);
445: PETSC_EXTERN PetscErrorCode DMPlexComputeRHSFunctionFVM(DM, PetscReal, Vec, Vec, void *);
446: PETSC_EXTERN PetscErrorCode DMPlexReconstructGradientsFVM(DM, Vec, Vec);
448: /* anchors */
449: PETSC_EXTERN PetscErrorCode DMPlexGetAnchors(DM, PetscSection *, IS *);
450: PETSC_EXTERN PetscErrorCode DMPlexSetAnchors(DM, PetscSection, IS);
451: /* tree */
452: PETSC_EXTERN PetscErrorCode DMPlexSetReferenceTree(DM, DM);
453: PETSC_EXTERN PetscErrorCode DMPlexGetReferenceTree(DM, DM *);
454: PETSC_EXTERN PetscErrorCode DMPlexReferenceTreeGetChildSymmetry(DM, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *);
455: PETSC_EXTERN PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm, PetscInt, PetscBool, DM *);
456: PETSC_EXTERN PetscErrorCode DMPlexSetTree(DM, PetscSection, PetscInt[], PetscInt[]);
457: PETSC_EXTERN PetscErrorCode DMPlexGetTree(DM, PetscSection *, PetscInt *[], PetscInt *[], PetscSection *, PetscInt *[]);
458: PETSC_EXTERN PetscErrorCode DMPlexGetTreeParent(DM, PetscInt, PetscInt *, PetscInt *);
459: PETSC_EXTERN PetscErrorCode DMPlexGetTreeChildren(DM, PetscInt, PetscInt *, const PetscInt *[]);
460: PETSC_EXTERN PetscErrorCode DMPlexTreeRefineCell(DM, PetscInt, DM *);
461: PETSC_EXTERN PetscErrorCode DMPlexComputeInjectorReferenceTree(DM, Mat *);
462: PETSC_EXTERN PetscErrorCode DMPlexTransferVecTree(DM, Vec, DM, Vec, PetscSF, PetscSF, PetscInt *, PetscInt *, PetscBool, PetscReal);
464: PETSC_EXTERN PetscErrorCode DMPlexMonitorThroughput(DM, void *);
466: /* natural order */
467: PETSC_EXTERN PetscErrorCode DMPlexCreateGlobalToNaturalSF(DM, PetscSection, PetscSF, PetscSF *);
468: PETSC_EXTERN PetscErrorCode DMPlexSetGlobalToNaturalSF(DM, PetscSF);
469: PETSC_EXTERN PetscErrorCode DMPlexGetGlobalToNaturalSF(DM, PetscSF *);
470: PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalBegin(DM, Vec, Vec);
471: PETSC_EXTERN PetscErrorCode DMPlexGlobalToNaturalEnd(DM, Vec, Vec);
472: PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalBegin(DM, Vec, Vec);
473: PETSC_EXTERN PetscErrorCode DMPlexNaturalToGlobalEnd(DM, Vec, Vec);
474: PETSC_EXTERN PetscErrorCode DMPlexCreateNaturalVector(DM, Vec *);
476: /* mesh adaptation */
477: PETSC_EXTERN PetscErrorCode DMPlexAdapt(DM, Vec, const char[], DM *);
478: PETSC_EXTERN PetscErrorCode DMPlexSnapToGeomModel(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
479: PETSC_EXTERN PetscErrorCode DMPlexMetricSetFromOptions(DM);
480: PETSC_EXTERN PetscErrorCode DMPlexMetricSetIsotropic(DM, PetscBool);
481: PETSC_EXTERN PetscErrorCode DMPlexMetricIsIsotropic(DM, PetscBool *);
482: PETSC_EXTERN PetscErrorCode DMPlexMetricSetUniform(DM, PetscBool);
483: PETSC_EXTERN PetscErrorCode DMPlexMetricIsUniform(DM, PetscBool *);
484: PETSC_EXTERN PetscErrorCode DMPlexMetricSetRestrictAnisotropyFirst(DM, PetscBool);
485: PETSC_EXTERN PetscErrorCode DMPlexMetricRestrictAnisotropyFirst(DM, PetscBool *);
486: PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoInsertion(DM, PetscBool);
487: PETSC_EXTERN PetscErrorCode DMPlexMetricNoInsertion(DM, PetscBool *);
488: PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoSwapping(DM, PetscBool);
489: PETSC_EXTERN PetscErrorCode DMPlexMetricNoSwapping(DM, PetscBool *);
490: PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoMovement(DM, PetscBool);
491: PETSC_EXTERN PetscErrorCode DMPlexMetricNoMovement(DM, PetscBool *);
492: PETSC_EXTERN PetscErrorCode DMPlexMetricSetNoSurf(DM, PetscBool);
493: PETSC_EXTERN PetscErrorCode DMPlexMetricNoSurf(DM, PetscBool *);
494: PETSC_EXTERN PetscErrorCode DMPlexMetricSetMinimumMagnitude(DM, PetscReal);
495: PETSC_EXTERN PetscErrorCode DMPlexMetricGetMinimumMagnitude(DM, PetscReal *);
496: PETSC_EXTERN PetscErrorCode DMPlexMetricSetMaximumMagnitude(DM, PetscReal);
497: PETSC_EXTERN PetscErrorCode DMPlexMetricGetMaximumMagnitude(DM, PetscReal *);
498: PETSC_EXTERN PetscErrorCode DMPlexMetricSetMaximumAnisotropy(DM, PetscReal);
499: PETSC_EXTERN PetscErrorCode DMPlexMetricGetMaximumAnisotropy(DM, PetscReal *);
500: PETSC_EXTERN PetscErrorCode DMPlexMetricSetTargetComplexity(DM, PetscReal);
501: PETSC_EXTERN PetscErrorCode DMPlexMetricGetTargetComplexity(DM, PetscReal *);
502: PETSC_EXTERN PetscErrorCode DMPlexMetricSetNormalizationOrder(DM, PetscReal);
503: PETSC_EXTERN PetscErrorCode DMPlexMetricGetNormalizationOrder(DM, PetscReal *);
504: PETSC_EXTERN PetscErrorCode DMPlexMetricSetGradationFactor(DM, PetscReal);
505: PETSC_EXTERN PetscErrorCode DMPlexMetricGetGradationFactor(DM, PetscReal *);
506: PETSC_EXTERN PetscErrorCode DMPlexMetricSetHausdorffNumber(DM, PetscReal);
507: PETSC_EXTERN PetscErrorCode DMPlexMetricGetHausdorffNumber(DM, PetscReal *);
508: PETSC_EXTERN PetscErrorCode DMPlexMetricSetVerbosity(DM, PetscInt);
509: PETSC_EXTERN PetscErrorCode DMPlexMetricGetVerbosity(DM, PetscInt *);
510: PETSC_EXTERN PetscErrorCode DMPlexMetricSetNumIterations(DM, PetscInt);
511: PETSC_EXTERN PetscErrorCode DMPlexMetricGetNumIterations(DM, PetscInt *);
512: PETSC_EXTERN PetscErrorCode DMPlexMetricCreate(DM, PetscInt, Vec *);
513: PETSC_EXTERN PetscErrorCode DMPlexMetricCreateUniform(DM, PetscInt, PetscReal, Vec *);
514: PETSC_EXTERN PetscErrorCode DMPlexMetricCreateIsotropic(DM, PetscInt, Vec, Vec *);
515: PETSC_EXTERN PetscErrorCode DMPlexMetricDeterminantCreate(DM, PetscInt, Vec *, DM *);
516: PETSC_EXTERN PetscErrorCode DMPlexMetricEnforceSPD(DM, Vec, PetscBool, PetscBool, Vec, Vec);
517: PETSC_EXTERN PetscErrorCode DMPlexMetricNormalize(DM, Vec, PetscBool, PetscBool, Vec, Vec);
518: PETSC_EXTERN PetscErrorCode DMPlexMetricAverage(DM, PetscInt, PetscReal[], Vec[], Vec);
519: PETSC_EXTERN PetscErrorCode DMPlexMetricAverage2(DM, Vec, Vec, Vec);
520: PETSC_EXTERN PetscErrorCode DMPlexMetricAverage3(DM, Vec, Vec, Vec, Vec);
521: PETSC_EXTERN PetscErrorCode DMPlexMetricIntersection(DM, PetscInt, Vec[], Vec);
522: PETSC_EXTERN PetscErrorCode DMPlexMetricIntersection2(DM, Vec, Vec, Vec);
523: PETSC_EXTERN PetscErrorCode DMPlexMetricIntersection3(DM, Vec, Vec, Vec, Vec);
525: PETSC_EXTERN PetscErrorCode DMPlexGlobalToLocalBasis(DM, Vec);
526: PETSC_EXTERN PetscErrorCode DMPlexLocalToGlobalBasis(DM, Vec);
527: PETSC_EXTERN PetscErrorCode DMPlexCreateBasisRotation(DM, PetscReal, PetscReal, PetscReal);
529: /* storage version */
530: #define DMPLEX_STORAGE_VERSION_FIRST "1.0.0"
531: #define DMPLEX_STORAGE_VERSION_STABLE "1.0.0"
532: #define DMPLEX_STORAGE_VERSION_LATEST "3.0.0"
534: PETSC_EXTERN PetscErrorCode DMPlexTopologyView(DM, PetscViewer);
535: PETSC_EXTERN PetscErrorCode DMPlexCoordinatesView(DM, PetscViewer);
536: PETSC_EXTERN PetscErrorCode DMPlexLabelsView(DM, PetscViewer);
537: PETSC_EXTERN PetscErrorCode DMPlexSectionView(DM, PetscViewer, DM);
538: PETSC_EXTERN PetscErrorCode DMPlexGlobalVectorView(DM, PetscViewer, DM, Vec);
539: PETSC_EXTERN PetscErrorCode DMPlexLocalVectorView(DM, PetscViewer, DM, Vec);
540: PETSC_EXTERN PetscErrorCode DMPlexVecView1D(DM, PetscInt, Vec[], PetscViewer);
542: PETSC_EXTERN PetscErrorCode DMPlexTopologyLoad(DM, PetscViewer, PetscSF *);
543: PETSC_EXTERN PetscErrorCode DMPlexCoordinatesLoad(DM, PetscViewer, PetscSF);
544: PETSC_EXTERN PetscErrorCode DMPlexLabelsLoad(DM, PetscViewer, PetscSF);
545: PETSC_EXTERN PetscErrorCode DMPlexSectionLoad(DM, PetscViewer, DM, PetscSF, PetscSF *, PetscSF *);
546: PETSC_EXTERN PetscErrorCode DMPlexGlobalVectorLoad(DM, PetscViewer, DM, PetscSF, Vec);
547: PETSC_EXTERN PetscErrorCode DMPlexLocalVectorLoad(DM, PetscViewer, DM, PetscSF, Vec);
549: PETSC_EXTERN PetscErrorCode DMPlexGetLocalOffsets(DM, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscInt **);
550: PETSC_EXTERN PetscErrorCode DMPlexGetLocalOffsetsSupport(DM, DMLabel, PetscInt, PetscInt *, PetscInt *, PetscInt *, PetscInt **, PetscInt **);
552: /* point queue */
553: PETSC_EXTERN PetscErrorCode DMPlexPointQueueCreate(PetscInt, DMPlexPointQueue *);
554: PETSC_EXTERN PetscErrorCode DMPlexPointQueueDestroy(DMPlexPointQueue *);
555: PETSC_EXTERN PetscErrorCode DMPlexPointQueueEnsureSize(DMPlexPointQueue);
556: PETSC_EXTERN PetscErrorCode DMPlexPointQueueEnqueue(DMPlexPointQueue, PetscInt);
557: PETSC_EXTERN PetscErrorCode DMPlexPointQueueDequeue(DMPlexPointQueue, PetscInt *);
558: PETSC_EXTERN PetscErrorCode DMPlexPointQueueFront(DMPlexPointQueue, PetscInt *);
559: PETSC_EXTERN PetscErrorCode DMPlexPointQueueBack(DMPlexPointQueue, PetscInt *);
560: PETSC_EXTERN PetscBool DMPlexPointQueueEmpty(DMPlexPointQueue);
561: PETSC_EXTERN PetscErrorCode DMPlexPointQueueEmptyCollective(PetscObject, DMPlexPointQueue, PetscBool *);
563: #endif