xref: /petsc/include/petscdmplex.h (revision 3e1910f1ab6113d8365e15c6b8c907ccce7ce4ea)
1 /*
2   DMPlex, for parallel unstructured distributed mesh problems.
3 */
4 #if !defined(__PETSCDMPLEX_H)
5 #define __PETSCDMPLEX_H
6 
7 #include <petscsftypes.h>
8 #include <petscdm.h>
9 
10 PETSC_EXTERN PetscErrorCode DMPlexCreate(MPI_Comm, DM*);
11 PETSC_EXTERN PetscErrorCode DMPlexCreateSubmesh(DM, const char[], DM*);
12 PETSC_EXTERN PetscErrorCode DMPlexCreateFromCellList(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscBool, const int[], PetscInt, const double[], DM*);
13 PETSC_EXTERN PetscErrorCode DMPlexCreateFromDAG(DM, PetscInt, const PetscInt [], const PetscInt [], const PetscInt [], const PetscInt [], const PetscScalar []);
14 PETSC_EXTERN PetscErrorCode DMPlexClone(DM, DM*);
15 PETSC_EXTERN PetscErrorCode DMPlexGetDimension(DM, PetscInt *);
16 PETSC_EXTERN PetscErrorCode DMPlexSetDimension(DM, PetscInt);
17 PETSC_EXTERN PetscErrorCode DMPlexGetChart(DM, PetscInt *, PetscInt *);
18 PETSC_EXTERN PetscErrorCode DMPlexSetChart(DM, PetscInt, PetscInt);
19 PETSC_EXTERN PetscErrorCode DMPlexGetConeSize(DM, PetscInt, PetscInt *);
20 PETSC_EXTERN PetscErrorCode DMPlexSetConeSize(DM, PetscInt, PetscInt);
21 PETSC_EXTERN PetscErrorCode DMPlexGetCone(DM, PetscInt, const PetscInt *[]);
22 PETSC_EXTERN PetscErrorCode DMPlexSetCone(DM, PetscInt, const PetscInt[]);
23 PETSC_EXTERN PetscErrorCode DMPlexInsertCone(DM, PetscInt, PetscInt, PetscInt);
24 PETSC_EXTERN PetscErrorCode DMPlexInsertConeOrientation(DM, PetscInt, PetscInt, PetscInt);
25 PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientation(DM, PetscInt, const PetscInt *[]);
26 PETSC_EXTERN PetscErrorCode DMPlexSetConeOrientation(DM, PetscInt, const PetscInt[]);
27 PETSC_EXTERN PetscErrorCode DMPlexGetSupportSize(DM, PetscInt, PetscInt *);
28 PETSC_EXTERN PetscErrorCode DMPlexSetSupportSize(DM, PetscInt, PetscInt);
29 PETSC_EXTERN PetscErrorCode DMPlexGetSupport(DM, PetscInt, const PetscInt *[]);
30 PETSC_EXTERN PetscErrorCode DMPlexSetSupport(DM, PetscInt, const PetscInt[]);
31 PETSC_EXTERN PetscErrorCode DMPlexInsertSupport(DM, PetscInt, PetscInt, PetscInt);
32 PETSC_EXTERN PetscErrorCode DMPlexGetConeSection(DM, PetscSection *);
33 PETSC_EXTERN PetscErrorCode DMPlexGetCones(DM, PetscInt *[]);
34 PETSC_EXTERN PetscErrorCode DMPlexGetConeOrientations(DM, PetscInt *[]);
35 PETSC_EXTERN PetscErrorCode DMPlexGetMaxSizes(DM, PetscInt *, PetscInt *);
36 PETSC_EXTERN PetscErrorCode DMPlexSymmetrize(DM);
37 PETSC_EXTERN PetscErrorCode DMPlexStratify(DM);
38 PETSC_EXTERN PetscErrorCode DMPlexGetCoordinateSection(DM, PetscSection *);
39 PETSC_EXTERN PetscErrorCode DMPlexSetCoordinateSection(DM, PetscSection);
40 PETSC_EXTERN PetscErrorCode DMPlexSetPreallocationCenterDimension(DM,PetscInt);
41 PETSC_EXTERN PetscErrorCode DMPlexPreallocateOperator(DM, PetscInt, PetscSection, PetscSection, PetscInt[], PetscInt[], PetscInt[], PetscInt[], Mat, PetscBool);
42 PETSC_EXTERN PetscErrorCode DMPlexGetPointLocal(DM,PetscInt,PetscInt*,PetscInt*);
43 PETSC_EXTERN PetscErrorCode DMPlexPointLocalRef(DM,PetscInt,PetscScalar*,void*);
44 PETSC_EXTERN PetscErrorCode DMPlexPointLocalRead(DM,PetscInt,const PetscScalar*,const void*);
45 PETSC_EXTERN PetscErrorCode DMPlexGetPointGlobal(DM,PetscInt,PetscInt*,PetscInt*);
46 PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRef(DM,PetscInt,PetscScalar*,void*);
47 PETSC_EXTERN PetscErrorCode DMPlexPointGlobalRead(DM,PetscInt,const PetscScalar*,const void*);
48 
49 /*S
50   DMLabel - Object which encapsulates a subset of the mesh from this DM
51 
52   Level: developer
53 
54   Concepts: grids, grid refinement
55 
56 .seealso:  DM, DMPlexCreate(), DMPlexCreateLabel()
57 S*/
58 typedef struct _n_DMLabel *DMLabel;
59 PETSC_EXTERN PetscErrorCode DMLabelCreate(const char [], DMLabel *);
60 PETSC_EXTERN PetscErrorCode DMLabelView(DMLabel, PetscViewer);
61 PETSC_EXTERN PetscErrorCode DMLabelDestroy(DMLabel *);
62 PETSC_EXTERN PetscErrorCode DMLabelDuplicate(DMLabel, DMLabel *);
63 PETSC_EXTERN PetscErrorCode DMLabelGetName(DMLabel, const char **);
64 PETSC_EXTERN PetscErrorCode DMLabelGetValue(DMLabel, PetscInt, PetscInt *);
65 PETSC_EXTERN PetscErrorCode DMLabelSetValue(DMLabel, PetscInt, PetscInt);
66 PETSC_EXTERN PetscErrorCode DMLabelClearValue(DMLabel, PetscInt, PetscInt);
67 PETSC_EXTERN PetscErrorCode DMLabelGetNumValues(DMLabel, PetscInt *);
68 PETSC_EXTERN PetscErrorCode DMLabelGetValueIS(DMLabel, IS *);
69 PETSC_EXTERN PetscErrorCode DMLabelGetStratumSize(DMLabel, PetscInt, PetscInt *);
70 PETSC_EXTERN PetscErrorCode DMLabelGetStratumIS(DMLabel, PetscInt, IS *);
71 PETSC_EXTERN PetscErrorCode DMLabelClearStratum(DMLabel, PetscInt);
72 PETSC_EXTERN PetscErrorCode DMLabelCreateIndex(DMLabel, PetscInt, PetscInt);
73 PETSC_EXTERN PetscErrorCode DMLabelDestroyIndex(DMLabel);
74 PETSC_EXTERN PetscErrorCode DMLabelHasPoint(DMLabel, PetscInt, PetscBool *);
75 
76 PETSC_EXTERN PetscErrorCode DMPlexCreateLabel(DM, const char []);
77 PETSC_EXTERN PetscErrorCode DMPlexGetLabelValue(DM, const char[], PetscInt, PetscInt *);
78 PETSC_EXTERN PetscErrorCode DMPlexSetLabelValue(DM, const char[], PetscInt, PetscInt);
79 PETSC_EXTERN PetscErrorCode DMPlexClearLabelValue(DM, const char[], PetscInt, PetscInt);
80 PETSC_EXTERN PetscErrorCode DMPlexGetLabelSize(DM, const char[], PetscInt *);
81 PETSC_EXTERN PetscErrorCode DMPlexGetLabelIdIS(DM, const char[], IS *);
82 PETSC_EXTERN PetscErrorCode DMPlexGetStratumSize(DM, const char [], PetscInt, PetscInt *);
83 PETSC_EXTERN PetscErrorCode DMPlexGetStratumIS(DM, const char [], PetscInt, IS *);
84 PETSC_EXTERN PetscErrorCode DMPlexClearLabelStratum(DM, const char[], PetscInt);
85 PETSC_EXTERN PetscErrorCode PetscSectionCreateGlobalSectionLabel(PetscSection, PetscSF, PetscBool, DMLabel, PetscInt, PetscSection *);
86 
87 PETSC_EXTERN PetscErrorCode DMPlexGetNumLabels(DM, PetscInt *);
88 PETSC_EXTERN PetscErrorCode DMPlexGetLabelName(DM, PetscInt, const char **);
89 PETSC_EXTERN PetscErrorCode DMPlexHasLabel(DM, const char [], PetscBool *);
90 PETSC_EXTERN PetscErrorCode DMPlexGetLabel(DM, const char *, DMLabel *);
91 PETSC_EXTERN PetscErrorCode DMPlexAddLabel(DM, DMLabel);
92 PETSC_EXTERN PetscErrorCode DMPlexRemoveLabel(DM, const char [], DMLabel *);
93 PETSC_EXTERN PetscErrorCode DMPlexGetCellNumbering(DM, IS *);
94 PETSC_EXTERN PetscErrorCode DMPlexGetVertexNumbering(DM, IS *);
95 
96 PETSC_EXTERN PetscErrorCode DMPlexGetMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
97 PETSC_EXTERN PetscErrorCode DMPlexGetFullMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
98 PETSC_EXTERN PetscErrorCode DMPlexRestoreMeet(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
99 PETSC_EXTERN PetscErrorCode DMPlexGetJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
100 PETSC_EXTERN PetscErrorCode DMPlexGetFullJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
101 PETSC_EXTERN PetscErrorCode DMPlexRestoreJoin(DM, PetscInt, const PetscInt [], PetscInt *, const PetscInt **);
102 PETSC_EXTERN PetscErrorCode DMPlexGetTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
103 PETSC_EXTERN PetscErrorCode DMPlexRestoreTransitiveClosure(DM, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
104 
105 PETSC_EXTERN PetscErrorCode DMPlexCreateNeighborCSR(DM, PetscInt, PetscInt *, PetscInt **, PetscInt **);
106 PETSC_EXTERN PetscErrorCode DMPlexCreatePartition(DM, PetscInt, PetscBool, PetscSection *, IS *, PetscSection *, IS *);
107 PETSC_EXTERN PetscErrorCode DMPlexCreatePartitionClosure(DM, PetscSection, IS, PetscSection *, IS *);
108 
109 PETSC_EXTERN PetscErrorCode DMPlexGenerate(DM, const char [], PetscBool , DM *);
110 PETSC_EXTERN PetscErrorCode DMPlexGetRefinementLimit(DM, PetscReal *);
111 PETSC_EXTERN PetscErrorCode DMPlexSetRefinementLimit(DM, PetscReal);
112 PETSC_EXTERN PetscErrorCode DMPlexGetRefinementUniform(DM, PetscBool *);
113 PETSC_EXTERN PetscErrorCode DMPlexSetRefinementUniform(DM, PetscBool);
114 PETSC_EXTERN PetscErrorCode DMPlexInterpolate(DM, DM *);
115 PETSC_EXTERN PetscErrorCode DMPlexDistribute(DM, const char[], PetscInt, DM*);
116 PETSC_EXTERN PetscErrorCode DMPlexLoad(PetscViewer, DM);
117 PETSC_EXTERN PetscErrorCode DMPlexGetSubpointMap(DM, DMLabel*);
118 PETSC_EXTERN PetscErrorCode DMPlexSetSubpointMap(DM, DMLabel);
119 PETSC_EXTERN PetscErrorCode DMPlexCreateSubpointIS(DM, IS *);
120 
121 PETSC_EXTERN PetscErrorCode DMPlexCreateCubeBoundary(DM, const PetscReal [], const PetscReal [], const PetscInt []);
122 PETSC_EXTERN PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm, PetscInt, PetscBool, DM *);
123 PETSC_EXTERN PetscErrorCode DMPlexCreateHexBoxMesh(MPI_Comm,PetscInt,const PetscInt[],DM *);
124 PETSC_EXTERN PetscErrorCode DMPlexGetDepth(DM, PetscInt *);
125 PETSC_EXTERN PetscErrorCode DMPlexGetDepthStratum(DM, PetscInt, PetscInt *, PetscInt *);
126 PETSC_EXTERN PetscErrorCode DMPlexGetHeightStratum(DM, PetscInt, PetscInt *, PetscInt *);
127 PETSC_EXTERN PetscErrorCode DMPlexCreateSection(DM, PetscInt, PetscInt,const PetscInt [],const PetscInt [], PetscInt,const PetscInt [],const IS [], PetscSection *);
128 PETSC_EXTERN PetscErrorCode DMPlexCreateConeSection(DM, PetscSection *);
129 PETSC_EXTERN PetscErrorCode DMPlexLabelCohesiveComplete(DM, DMLabel);
130 PETSC_EXTERN PetscErrorCode DMPlexMarkBoundaryFaces(DM, DMLabel);
131 PETSC_EXTERN PetscErrorCode DMPlexLabelComplete(DM, DMLabel);
132 
133 /* Support for cell-vertex meshes */
134 PETSC_EXTERN PetscErrorCode DMPlexGetNumFaceVertices(DM, PetscInt, PetscInt, PetscInt *);
135 
136 /* FEM Support */
137 PETSC_EXTERN PetscErrorCode DMPlexComputeCellGeometry(DM, PetscInt, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
138 PETSC_EXTERN PetscErrorCode DMPlexVecGetClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]);
139 PETSC_EXTERN PetscErrorCode DMPlexVecRestoreClosure(DM, PetscSection, Vec, PetscInt, PetscInt *, PetscScalar *[]);
140 PETSC_EXTERN PetscErrorCode DMPlexVecSetClosure(DM, PetscSection, Vec, PetscInt, const PetscScalar[], InsertMode);
141 PETSC_EXTERN PetscErrorCode DMPlexMatSetClosure(DM, PetscSection, PetscSection, Mat, PetscInt, const PetscScalar[], InsertMode);
142 
143 PETSC_EXTERN PetscErrorCode DMPlexCreateExodus(MPI_Comm, PetscInt, PetscBool, DM *);
144 
145 PETSC_EXTERN PetscErrorCode DMPlexConstructGhostCells(DM, const char [], PetscInt *, DM *);
146 PETSC_EXTERN PetscErrorCode DMPlexConstructCohesiveCells(DM, DMLabel, DM *);
147 
148 PETSC_EXTERN PetscErrorCode DMPlexGetHybridBounds(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
149 PETSC_EXTERN PetscErrorCode DMPlexSetHybridBounds(DM, PetscInt, PetscInt, PetscInt, PetscInt);
150 PETSC_EXTERN PetscErrorCode DMPlexGetVTKCellHeight(DM, PetscInt *);
151 PETSC_EXTERN PetscErrorCode DMPlexSetVTKCellHeight(DM, PetscInt);
152 PETSC_EXTERN PetscErrorCode DMPlexVTKWriteAll(PetscObject, PetscViewer);
153 
154 PETSC_EXTERN PetscErrorCode DMPlexGetScale(DM, PetscUnit, PetscReal *);
155 PETSC_EXTERN PetscErrorCode DMPlexSetScale(DM, PetscUnit, PetscReal);
156 
157 typedef struct {
158   DM    dm;
159   Vec   u; /* The base vector for the Jacbobian action J(u) x */
160   Mat   J; /* Preconditioner for testing */
161   void *user;
162 } JacActionCtx;
163 
164 PETSC_EXTERN PetscErrorCode DMPlexProjectFunction(DM, PetscInt, PetscScalar (**)(const PetscReal []), InsertMode, Vec);
165 PETSC_EXTERN PetscErrorCode DMPlexProjectFunctionLocal(DM, PetscInt, PetscScalar (**)(const PetscReal []), InsertMode, Vec);
166 PETSC_EXTERN PetscErrorCode DMPlexComputeL2Diff(DM, PetscQuadrature[], PetscScalar (**)(const PetscReal []), Vec, PetscReal *);
167 PETSC_EXTERN PetscErrorCode DMPlexComputeResidualFEM(DM, Vec, Vec, void *);
168 PETSC_EXTERN PetscErrorCode DMPlexComputeJacobianActionFEM(DM, Mat, Vec, Vec, void *);
169 PETSC_EXTERN PetscErrorCode DMPlexComputeJacobianFEM(DM, Vec, Mat, Mat, MatStructure*,void *);
170 PETSC_EXTERN PetscErrorCode DMPlexSetFEMIntegration(DM,
171                                                        PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[],
172                                                                           const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[],
173                                                                           void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
174                                                                           void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]),
175                                                        PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], const PetscScalar[],
176                                                                           const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[],
177                                                                           void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
178                                                                           void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
179                                                                           void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
180                                                                           void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]),
181                                                        PetscErrorCode (*)(PetscInt, PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[],
182                                                                           const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[],
183                                                                           void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
184                                                                           void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
185                                                                           void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
186                                                                           void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]));
187 #endif
188