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