xref: /petsc/include/petscdm.h (revision 0b4b7b1c20c2ed4ade67e3d50a7710fe0ffbfca5)
1 /*
2       Objects to manage the interactions between the mesh data structures and the algebraic objects
3 */
4 #pragma once
5 #include "petscsystypes.h"
6 #include <petscmat.h>
7 #include <petscdmtypes.h>
8 #include <petscfetypes.h>
9 #include <petscdstypes.h>
10 #include <petscdmlabel.h>
11 #include <petscdt.h>
12 
13 /* SUBMANSEC = DM */
14 
15 PETSC_EXTERN PetscErrorCode DMInitializePackage(void);
16 
17 PETSC_EXTERN PetscClassId DM_CLASSID;
18 
19 #define DMLOCATEPOINT_POINT_NOT_FOUND -367
20 
21 /*J
22    DMType - String with the name of a PETSc `DM`. These are all the `DM` provided by PETSc.
23 
24    Level: beginner
25 
26    Note:
27    These can be used with `DMSetType()` or the options database key `-dm_type` to set the specific data structures and algorithms to use with a specific `DM`.
28    But more commonly one calls directly a constructor for a particular `DMType` such as `DMDACreate()`
29 
30 .seealso: [](ch_dmbase), `DMSetType()`, `DMCreate()`, `DM`, `DMDACreate()`
31 J*/
32 typedef const char *DMType;
33 #define DMDA        "da"
34 #define DMCOMPOSITE "composite"
35 #define DMSLICED    "sliced"
36 #define DMSHELL     "shell"
37 #define DMPLEX      "plex"
38 #define DMREDUNDANT "redundant"
39 #define DMPATCH     "patch"
40 #define DMMOAB      "moab"
41 #define DMNETWORK   "network"
42 #define DMFOREST    "forest"
43 #define DMP4EST     "p4est"
44 #define DMP8EST     "p8est"
45 #define DMSWARM     "swarm"
46 #define DMPRODUCT   "product"
47 #define DMSTAG      "stag"
48 
49 PETSC_EXTERN const char *const       DMBoundaryTypes[];
50 PETSC_EXTERN const char *const       DMBoundaryConditionTypes[];
51 PETSC_EXTERN const char *const       DMBlockingTypes[];
52 PETSC_EXTERN PetscFunctionList       DMList;
53 PETSC_EXTERN DMGeneratorFunctionList DMGenerateList;
54 PETSC_EXTERN PetscFunctionList       DMGeomModelList;
55 PETSC_EXTERN PetscErrorCode          DMCreate(MPI_Comm, DM *);
56 PETSC_EXTERN PetscErrorCode          DMClone(DM, DM *);
57 PETSC_EXTERN PetscErrorCode          DMSetType(DM, DMType);
58 PETSC_EXTERN PetscErrorCode          DMGetType(DM, DMType *);
59 PETSC_EXTERN PetscErrorCode          DMRegister(const char[], PetscErrorCode (*)(DM));
60 PETSC_EXTERN PetscErrorCode          DMRegisterDestroy(void);
61 
62 PETSC_EXTERN PetscErrorCode DMView(DM, PetscViewer);
63 PETSC_EXTERN PetscErrorCode DMLoad(DM, PetscViewer);
64 PETSC_EXTERN PetscErrorCode DMDestroy(DM *);
65 PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM, Vec *);
66 PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM, Vec *);
67 PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM, Vec *);
68 PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM, Vec *);
69 PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM, Vec *);
70 PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM, Vec *);
71 PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
72 PETSC_EXTERN PetscErrorCode DMClearLocalVectors(DM);
73 PETSC_EXTERN PetscErrorCode DMClearNamedGlobalVectors(DM);
74 PETSC_EXTERN PetscErrorCode DMClearNamedLocalVectors(DM);
75 PETSC_EXTERN PetscErrorCode DMHasNamedGlobalVector(DM, const char *, PetscBool *);
76 PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM, const char *, Vec *);
77 PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM, const char *, Vec *);
78 PETSC_EXTERN PetscErrorCode DMHasNamedLocalVector(DM, const char *, PetscBool *);
79 PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM, const char *, Vec *);
80 PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM, const char *, Vec *);
81 PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM, ISLocalToGlobalMapping *);
82 PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM, PetscInt *, char ***, IS **);
83 PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM, PetscInt *);
84 PETSC_EXTERN PetscErrorCode DMCreateColoring(DM, ISColoringType, ISColoring *);
85 PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM, Mat *);
86 PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateSkip(DM, PetscBool);
87 PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM, PetscBool);
88 PETSC_EXTERN PetscErrorCode DMSetMatrixStructureOnly(DM, PetscBool);
89 PETSC_EXTERN PetscErrorCode DMSetBlockingType(DM, DMBlockingType);
90 PETSC_EXTERN PetscErrorCode DMGetBlockingType(DM, DMBlockingType *);
91 PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM, DM, Mat *, Vec *);
92 PETSC_EXTERN PetscErrorCode DMCreateRestriction(DM, DM, Mat *);
93 PETSC_EXTERN PetscErrorCode DMCreateInjection(DM, DM, Mat *);
94 PETSC_EXTERN PetscErrorCode DMCreateMassMatrix(DM, DM, Mat *);
95 PETSC_EXTERN PetscErrorCode DMCreateMassMatrixLumped(DM, Vec *, Vec *);
96 PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM, PetscInt, MPI_Datatype, void *);
97 PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM, PetscInt, MPI_Datatype, void *);
98 PETSC_EXTERN PetscErrorCode DMRefine(DM, MPI_Comm, DM *);
99 PETSC_EXTERN PetscErrorCode DMCoarsen(DM, MPI_Comm, DM *);
100 PETSC_EXTERN PetscErrorCode DMGetCoarseDM(DM, DM *);
101 PETSC_EXTERN PetscErrorCode DMSetCoarseDM(DM, DM);
102 PETSC_EXTERN PetscErrorCode DMGetFineDM(DM, DM *);
103 PETSC_EXTERN PetscErrorCode DMSetFineDM(DM, DM);
104 PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM, PetscInt, DM[]);
105 PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM, PetscInt, DM[]);
106 PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, Vec, Mat, DM, void *), void *);
107 PETSC_EXTERN PetscErrorCode DMCoarsenHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, Vec, Mat, DM, void *), void *);
108 PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, DM, void *), void *);
109 PETSC_EXTERN PetscErrorCode DMRefineHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, DM, void *), void *);
110 PETSC_EXTERN PetscErrorCode DMRestrict(DM, Mat, Vec, Mat, DM);
111 PETSC_EXTERN PetscErrorCode DMInterpolate(DM, Mat, DM);
112 PETSC_EXTERN PetscErrorCode DMInterpolateSolution(DM, DM, Mat, Vec, Vec);
113 PETSC_EXTERN PetscErrorCode DMExtrude(DM, PetscInt, DM *);
114 PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
115 PETSC_EXTERN PetscErrorCode DMViewFromOptions(DM, PetscObject, const char[]);
116 
117 PETSC_EXTERN PetscErrorCode DMGenerate(DM, const char[], PetscBool, DM *);
118 PETSC_EXTERN PetscErrorCode DMGenerateRegister(const char[], PetscErrorCode (*)(DM, PetscBool, DM *), PetscErrorCode (*)(DM, PetscReal *, DM *), PetscErrorCode (*)(DM, Vec, DMLabel, DMLabel, DM *), PetscInt);
119 PETSC_EXTERN PetscErrorCode DMGenerateRegisterAll(void);
120 PETSC_EXTERN PetscErrorCode DMGenerateRegisterDestroy(void);
121 PETSC_EXTERN PetscErrorCode DMGeomModelRegister(const char[], PetscErrorCode (*)(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]));
122 PETSC_EXTERN PetscErrorCode DMGeomModelRegisterAll(void);
123 PETSC_EXTERN PetscErrorCode DMGeomModelRegisterDestroy(void);
124 PETSC_EXTERN PetscErrorCode DMAdaptLabel(DM, DMLabel, DM *);
125 PETSC_EXTERN PetscErrorCode DMAdaptMetric(DM, Vec, DMLabel, DMLabel, DM *);
126 
127 PETSC_EXTERN PetscErrorCode DMSetUp(DM);
128 PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM, DM, Mat, Vec *);
129 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMDACreateAggregates()", ) PetscErrorCode DMCreateAggregates(DM, DM, Mat *);
130 PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM, PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), void *);
131 PETSC_EXTERN PetscErrorCode DMLocalToGlobalHookAdd(DM, PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), void *);
132 PETSC_EXTERN PetscErrorCode DMGlobalToLocal(DM, Vec, InsertMode, Vec);
133 PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM, Vec, InsertMode, Vec);
134 PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM, Vec, InsertMode, Vec);
135 PETSC_EXTERN PetscErrorCode DMLocalToGlobal(DM, Vec, InsertMode, Vec);
136 PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM, Vec, InsertMode, Vec);
137 PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM, Vec, InsertMode, Vec);
138 PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM, Vec, InsertMode, Vec);
139 PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM, Vec, InsertMode, Vec);
140 PETSC_EXTERN PetscErrorCode DMConvert(DM, DMType, DM *);
141 
142 /* Topology support */
143 PETSC_EXTERN PetscErrorCode DMGetDimension(DM, PetscInt *);
144 PETSC_EXTERN PetscErrorCode DMSetDimension(DM, PetscInt);
145 PETSC_EXTERN PetscErrorCode DMGetDimPoints(DM, PetscInt, PetscInt *, PetscInt *);
146 PETSC_EXTERN PetscErrorCode DMGetUseNatural(DM, PetscBool *);
147 PETSC_EXTERN PetscErrorCode DMSetUseNatural(DM, PetscBool);
148 PETSC_EXTERN PetscErrorCode DMGetNeighbors(DM, PetscInt *, const PetscMPIInt **);
149 
150 /* Coordinate support */
151 PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM, DM *);
152 PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM, DM);
153 PETSC_EXTERN PetscErrorCode DMGetCellCoordinateDM(DM, DM *);
154 PETSC_EXTERN PetscErrorCode DMSetCellCoordinateDM(DM, DM);
155 PETSC_EXTERN PetscErrorCode DMGetCoordinateDim(DM, PetscInt *);
156 PETSC_EXTERN PetscErrorCode DMSetCoordinateDim(DM, PetscInt);
157 PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM, PetscSection *);
158 PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM, PetscInt, PetscSection);
159 PETSC_EXTERN PetscErrorCode DMGetCellCoordinateSection(DM, PetscSection *);
160 PETSC_EXTERN PetscErrorCode DMSetCellCoordinateSection(DM, PetscInt, PetscSection);
161 PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM, Vec *);
162 PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM, Vec);
163 PETSC_EXTERN PetscErrorCode DMGetCellCoordinates(DM, Vec *);
164 PETSC_EXTERN PetscErrorCode DMSetCellCoordinates(DM, Vec);
165 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalSetUp(DM);
166 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM, Vec *);
167 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalNoncollective(DM, Vec *);
168 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalTuple(DM, IS, PetscSection *, Vec *);
169 PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM, Vec);
170 PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM);
171 PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocal(DM, Vec *);
172 PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM, Vec *);
173 PETSC_EXTERN PetscErrorCode DMSetCellCoordinatesLocal(DM, Vec);
174 PETSC_EXTERN PetscErrorCode DMGetCoordinateField(DM, DMField *);
175 PETSC_EXTERN PetscErrorCode DMSetCoordinateField(DM, DMField);
176 PETSC_EXTERN PetscErrorCode DMGetLocalBoundingBox(DM, PetscReal[], PetscReal[]);
177 PETSC_EXTERN PetscErrorCode DMGetBoundingBox(DM, PetscReal[], PetscReal[]);
178 PETSC_EXTERN PetscErrorCode DMSetCoordinateDisc(DM, PetscFE, PetscBool);
179 PETSC_EXTERN PetscErrorCode DMLocatePoints(DM, Vec, DMPointLocationType, PetscSF *);
180 PETSC_EXTERN PetscErrorCode DMSnapToGeomModel(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
181 PETSC_EXTERN PetscErrorCode DMSetSnapToGeomModel(DM, const char[]);
182 
183 /* Periodicity support */
184 PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM, const PetscReal *[], const PetscReal *[], const PetscReal *[]);
185 PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM, const PetscReal[], const PetscReal[], const PetscReal[]);
186 PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate(DM, const PetscScalar[], PetscBool, PetscScalar[]);
187 PETSC_EXTERN PetscErrorCode DMLocalizeCoordinates(DM);
188 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalized(DM, PetscBool *);
189 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalizedLocal(DM, PetscBool *);
190 PETSC_EXTERN PetscErrorCode DMGetSparseLocalize(DM, PetscBool *);
191 PETSC_EXTERN PetscErrorCode DMSetSparseLocalize(DM, PetscBool);
192 
193 /* block hook interface */
194 PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, VecScatter, VecScatter, DM, void *), void *);
195 PETSC_EXTERN PetscErrorCode DMSubDomainHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, VecScatter, VecScatter, DM, void *), void *);
196 PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM, VecScatter, VecScatter, DM);
197 
198 PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM, const char[]);
199 PETSC_EXTERN PetscErrorCode DMAppendOptionsPrefix(DM, const char[]);
200 PETSC_EXTERN PetscErrorCode DMGetOptionsPrefix(DM, const char *[]);
201 PETSC_EXTERN PetscErrorCode DMSetVecType(DM, VecType);
202 PETSC_EXTERN PetscErrorCode DMGetVecType(DM, VecType *);
203 PETSC_EXTERN PetscErrorCode DMSetMatType(DM, MatType);
204 PETSC_EXTERN PetscErrorCode DMGetMatType(DM, MatType *);
205 PETSC_EXTERN PetscErrorCode DMSetISColoringType(DM, ISColoringType);
206 PETSC_EXTERN PetscErrorCode DMGetISColoringType(DM, ISColoringType *);
207 PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM, void *);
208 PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM, PetscErrorCode (*)(void **));
209 PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM, void *);
210 PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM, PetscErrorCode (*)(DM, Vec, Vec));
211 PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM, PetscBool *);
212 PETSC_EXTERN PetscErrorCode DMHasColoring(DM, PetscBool *);
213 PETSC_EXTERN PetscErrorCode DMHasCreateRestriction(DM, PetscBool *);
214 PETSC_EXTERN PetscErrorCode DMHasCreateInjection(DM, PetscBool *);
215 PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM, Vec, Vec);
216 
217 PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, const PetscInt[], IS *, DM *);
218 PETSC_EXTERN PetscErrorCode DMCreateSuperDM(DM[], PetscInt, IS **, DM *);
219 PETSC_EXTERN PetscErrorCode DMCreateSectionSubDM(DM, PetscInt, const PetscInt[], const PetscInt[], const PetscInt[], IS *, DM *);
220 PETSC_EXTERN PetscErrorCode DMCreateSectionSuperDM(DM[], PetscInt, IS **, DM *);
221 PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM, PetscInt *, char ***, IS **, DM **);
222 PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM, PetscInt *, char ***, IS **, IS **, DM **);
223 PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM, PetscInt, DM *, VecScatter **, VecScatter **, VecScatter **);
224 
225 PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM, PetscInt *);
226 PETSC_EXTERN PetscErrorCode DMSetRefineLevel(DM, PetscInt);
227 PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM, PetscInt *);
228 PETSC_EXTERN PetscErrorCode DMSetCoarsenLevel(DM, PetscInt);
229 PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);
230 
231 PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM *);
232 PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
233 PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM *);
234 PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
235 PETSC_EXTERN PetscErrorCode MatFDColoringUseDM(Mat, MatFDColoring);
236 
237 typedef struct NLF_DAAD *NLF;
238 
239 #define DM_FILE_CLASSID 1211221
240 
241 /* FEM support */
242 PETSC_EXTERN PetscErrorCode DMPrintCellIndices(PetscInt, const char[], PetscInt, const PetscInt[]);
243 PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char[], PetscInt, const PetscScalar[]);
244 PETSC_EXTERN PetscErrorCode DMPrintCellVectorReal(PetscInt, const char[], PetscInt, const PetscReal[]);
245 PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char[], PetscInt, PetscInt, const PetscScalar[]);
246 PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char[], PetscReal, Vec);
247 
248 PETSC_EXTERN PetscErrorCode DMSetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *));
249 PETSC_EXTERN PetscErrorCode DMGetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *));
250 PETSC_EXTERN PetscErrorCode DMSetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *));
251 PETSC_EXTERN PetscErrorCode DMGetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *));
252 
253 PETSC_EXTERN PetscErrorCode DMGetSection(DM, PetscSection *); /* Use DMGetLocalSection() in new code (since v3.12) */
254 PETSC_EXTERN PetscErrorCode DMSetSection(DM, PetscSection);   /* Use DMSetLocalSection() in new code (since v3.12) */
255 PETSC_EXTERN PetscErrorCode DMGetLocalSection(DM, PetscSection *);
256 PETSC_EXTERN PetscErrorCode DMSetLocalSection(DM, PetscSection);
257 PETSC_EXTERN PetscErrorCode DMGetGlobalSection(DM, PetscSection *);
258 PETSC_EXTERN PetscErrorCode DMSetGlobalSection(DM, PetscSection);
259 PETSC_EXTERN PetscErrorCode DMCreateSectionPermutation(DM, IS *, PetscBT *);
260 PETSC_EXTERN PetscErrorCode DMReorderSectionGetDefault(DM, DMReorderDefaultFlag *);
261 PETSC_EXTERN PetscErrorCode DMReorderSectionSetDefault(DM, DMReorderDefaultFlag);
262 PETSC_EXTERN PetscErrorCode DMReorderSectionGetType(DM, MatOrderingType *);
263 PETSC_EXTERN PetscErrorCode DMReorderSectionSetType(DM, MatOrderingType);
264 PETSC_EXTERN PetscErrorCode DMUseTensorOrder(DM, PetscBool);
265 static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMGetSection()", ) PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *s)
266 {
267   return DMGetSection(dm, s);
268 }
269 static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMSetSection()", ) PetscErrorCode DMSetDefaultSection(DM dm, PetscSection s)
270 {
271   return DMSetSection(dm, s);
272 }
273 static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMGetGlobalSection()", ) PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *s)
274 {
275   return DMGetGlobalSection(dm, s);
276 }
277 static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMSetGlobalSection()", ) PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection s)
278 {
279   return DMSetGlobalSection(dm, s);
280 }
281 
282 PETSC_EXTERN PetscErrorCode DMGetSectionSF(DM, PetscSF *);
283 PETSC_EXTERN PetscErrorCode DMSetSectionSF(DM, PetscSF);
284 PETSC_EXTERN PetscErrorCode DMCreateSectionSF(DM, PetscSection, PetscSection);
285 static inline PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMGetSectionSF()", ) PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *s)
286 {
287   return DMGetSectionSF(dm, s);
288 }
289 static inline PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMSetSectionSF()", ) PetscErrorCode DMSetDefaultSF(DM dm, PetscSF s)
290 {
291   return DMSetSectionSF(dm, s);
292 }
293 static inline PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMCreateSectionSF()", ) PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection l, PetscSection g)
294 {
295   return DMCreateSectionSF(dm, l, g);
296 }
297 PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
298 PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);
299 PETSC_EXTERN PetscErrorCode DMGetNaturalSF(DM, PetscSF *);
300 PETSC_EXTERN PetscErrorCode DMSetNaturalSF(DM, PetscSF);
301 
302 PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *, Vec *);
303 PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat, Vec);
304 
305 PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *);
306 PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *);
307 PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal);
308 PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char[], PetscInt, PetscReal *);
309 PETSC_EXTERN PetscErrorCode DMGetOutputSequenceLength(DM, PetscViewer, const char[], PetscInt *);
310 
311 PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
312 PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
313 PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, DMLabel *, PetscObject *);
314 PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, DMLabel, PetscObject);
315 PETSC_EXTERN PetscErrorCode DMAddField(DM, DMLabel, PetscObject);
316 PETSC_EXTERN PetscErrorCode DMSetFieldAvoidTensor(DM, PetscInt, PetscBool);
317 PETSC_EXTERN PetscErrorCode DMGetFieldAvoidTensor(DM, PetscInt, PetscBool *);
318 PETSC_EXTERN PetscErrorCode DMClearFields(DM);
319 PETSC_EXTERN PetscErrorCode DMCopyFields(DM, PetscInt, PetscInt, DM);
320 PETSC_EXTERN PetscErrorCode DMGetAdjacency(DM, PetscInt, PetscBool *, PetscBool *);
321 PETSC_EXTERN PetscErrorCode DMSetAdjacency(DM, PetscInt, PetscBool, PetscBool);
322 PETSC_EXTERN PetscErrorCode DMGetBasicAdjacency(DM, PetscBool *, PetscBool *);
323 PETSC_EXTERN PetscErrorCode DMSetBasicAdjacency(DM, PetscBool, PetscBool);
324 
325 PETSC_EXTERN PetscErrorCode DMGetNumDS(DM, PetscInt *);
326 PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *);
327 PETSC_EXTERN PetscErrorCode DMGetCellDS(DM, PetscInt, PetscDS *, PetscDS *);
328 PETSC_EXTERN PetscErrorCode DMGetRegionDS(DM, DMLabel, IS *, PetscDS *, PetscDS *);
329 PETSC_EXTERN PetscErrorCode DMSetRegionDS(DM, DMLabel, IS, PetscDS, PetscDS);
330 PETSC_EXTERN PetscErrorCode DMGetRegionNumDS(DM, PetscInt, DMLabel *, IS *, PetscDS *, PetscDS *);
331 PETSC_EXTERN PetscErrorCode DMSetRegionNumDS(DM, PetscInt, DMLabel, IS, PetscDS, PetscDS);
332 PETSC_EXTERN PetscErrorCode DMFindRegionNum(DM, PetscDS, PetscInt *);
333 PETSC_EXTERN PetscErrorCode DMCreateFEDefault(DM, PetscInt, const char[], PetscInt, PetscFE *);
334 PETSC_EXTERN PetscErrorCode DMCreateDS(DM);
335 PETSC_EXTERN PetscErrorCode DMClearDS(DM);
336 PETSC_EXTERN PetscErrorCode DMCopyDS(DM, PetscInt, PetscInt, DM);
337 PETSC_EXTERN PetscErrorCode DMCopyDisc(DM, DM);
338 PETSC_EXTERN PetscErrorCode DMComputeExactSolution(DM, PetscReal, Vec, Vec);
339 PETSC_EXTERN PetscErrorCode DMGetNumAuxiliaryVec(DM, PetscInt *);
340 PETSC_EXTERN PetscErrorCode DMGetAuxiliaryVec(DM, DMLabel, PetscInt, PetscInt, Vec *);
341 PETSC_EXTERN PetscErrorCode DMSetAuxiliaryVec(DM, DMLabel, PetscInt, PetscInt, Vec);
342 PETSC_EXTERN PetscErrorCode DMGetAuxiliaryLabels(DM, DMLabel[], PetscInt[], PetscInt[]);
343 PETSC_EXTERN PetscErrorCode DMCopyAuxiliaryVec(DM, DM);
344 PETSC_EXTERN PetscErrorCode DMClearAuxiliaryVec(DM);
345 
346 /*MC
347   DMInterpolationInfo - Structure for holding information about interpolation on a mesh
348 
349   Synopsis:
350     comm   - The communicator
351     dim    - The spatial dimension of points
352     nInput - The number of input points
353     points - The input point coordinates
354     cells  - The cell containing each point
355     n      - The number of local points
356     coords - The point coordinates
357     dof    - The number of components to interpolate
358 
359   Level: intermediate
360 
361 .seealso: [](ch_dmbase), `DM`, `DMInterpolationCreate()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
362 M*/
363 struct _DMInterpolationInfo {
364   MPI_Comm   comm;
365   PetscInt   dim;    /*1 The spatial dimension of points */
366   PetscInt   nInput; /* The number of input points */
367   PetscReal *points; /* The input point coordinates */
368   PetscInt  *cells;  /* The cell containing each point */
369   PetscInt   n;      /* The number of local points */
370   Vec        coords; /* The point coordinates */
371   PetscInt   dof;    /* The number of components to interpolate */
372 };
373 typedef struct _DMInterpolationInfo *DMInterpolationInfo;
374 
375 PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
376 PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
377 PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
378 PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
379 PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
380 PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
381 PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool, PetscBool);
382 PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
383 PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
384 PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
385 PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
386 PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);
387 
388 PETSC_EXTERN PetscErrorCode DMCreateLabel(DM, const char[]);
389 PETSC_EXTERN PetscErrorCode DMCreateLabelAtIndex(DM, PetscInt, const char[]);
390 PETSC_EXTERN PetscErrorCode DMGetLabelValue(DM, const char[], PetscInt, PetscInt *);
391 PETSC_EXTERN PetscErrorCode DMSetLabelValue(DM, const char[], PetscInt, PetscInt);
392 PETSC_EXTERN PetscErrorCode DMClearLabelValue(DM, const char[], PetscInt, PetscInt);
393 PETSC_EXTERN PetscErrorCode DMGetLabelSize(DM, const char[], PetscInt *);
394 PETSC_EXTERN PetscErrorCode DMGetLabelIdIS(DM, const char[], IS *);
395 PETSC_EXTERN PetscErrorCode DMGetStratumSize(DM, const char[], PetscInt, PetscInt *);
396 PETSC_EXTERN PetscErrorCode DMGetStratumIS(DM, const char[], PetscInt, IS *);
397 PETSC_EXTERN PetscErrorCode DMSetStratumIS(DM, const char[], PetscInt, IS);
398 PETSC_EXTERN PetscErrorCode DMClearLabelStratum(DM, const char[], PetscInt);
399 PETSC_EXTERN PetscErrorCode DMGetLabelOutput(DM, const char[], PetscBool *);
400 PETSC_EXTERN PetscErrorCode DMSetLabelOutput(DM, const char[], PetscBool);
401 PETSC_EXTERN PetscErrorCode DMGetFirstLabeledPoint(DM, DM, DMLabel, PetscInt, const PetscInt *, PetscInt, PetscInt *, PetscDS *);
402 
403 /*E
404    DMCopyLabelsMode - Determines how `DMCopyLabels()` behaves when there is a `DMLabel` in the source and destination `DM`s with the same name
405 
406    Values:
407 +  `DM_COPY_LABELS_REPLACE` - replace label in destination by label from source
408 .  `DM_COPY_LABELS_KEEP`    - keep destination label
409 -  `DM_COPY_LABELS_FAIL`    - generate an error
410 
411    Level: advanced
412 
413 .seealso: [](ch_dmbase), `DMLabel`, `DM`, `DMCompareLabels()`, `DMRemoveLabel()`
414 E*/
415 typedef enum {
416   DM_COPY_LABELS_REPLACE,
417   DM_COPY_LABELS_KEEP,
418   DM_COPY_LABELS_FAIL
419 } DMCopyLabelsMode;
420 PETSC_EXTERN const char *const DMCopyLabelsModes[];
421 
422 PETSC_EXTERN PetscErrorCode DMGetNumLabels(DM, PetscInt *);
423 PETSC_EXTERN PetscErrorCode DMGetLabelName(DM, PetscInt, const char **);
424 PETSC_EXTERN PetscErrorCode DMHasLabel(DM, const char[], PetscBool *);
425 PETSC_EXTERN PetscErrorCode DMGetLabel(DM, const char *, DMLabel *);
426 PETSC_EXTERN PetscErrorCode DMSetLabel(DM, DMLabel);
427 PETSC_EXTERN PetscErrorCode DMGetLabelByNum(DM, PetscInt, DMLabel *);
428 PETSC_EXTERN PetscErrorCode DMAddLabel(DM, DMLabel);
429 PETSC_EXTERN PetscErrorCode DMRemoveLabel(DM, const char[], DMLabel *);
430 PETSC_EXTERN PetscErrorCode DMRemoveLabelBySelf(DM, DMLabel *, PetscBool);
431 PETSC_EXTERN PetscErrorCode DMCopyLabels(DM, DM, PetscCopyMode, PetscBool, DMCopyLabelsMode emode);
432 PETSC_EXTERN PetscErrorCode DMCompareLabels(DM, DM, PetscBool *, char **);
433 
434 PETSC_EXTERN PetscErrorCode DMAddBoundary(DM, DMBoundaryConditionType, const char[], DMLabel, PetscInt, const PetscInt[], PetscInt, PetscInt, const PetscInt[], void (*)(void), void (*)(void), void *, PetscInt *);
435 PETSC_EXTERN PetscErrorCode DMIsBoundaryPoint(DM, PetscInt, PetscBool *);
436 
437 PETSC_EXTERN PetscErrorCode DMProjectFunction(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
438 PETSC_EXTERN PetscErrorCode DMProjectFunctionLocal(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
439 PETSC_EXTERN PetscErrorCode DMProjectFunctionLabel(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
440 PETSC_EXTERN PetscErrorCode DMProjectFunctionLabelLocal(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
441 PETSC_EXTERN PetscErrorCode DMProjectFieldLocal(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);
442 PETSC_EXTERN PetscErrorCode DMProjectFieldLabel(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Vec, void (**funcs)(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);
443 PETSC_EXTERN PetscErrorCode DMProjectFieldLabelLocal(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);
444 PETSC_EXTERN PetscErrorCode DMProjectBdFieldLabelLocal(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);
445 PETSC_EXTERN PetscErrorCode DMComputeL2Diff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
446 PETSC_EXTERN PetscErrorCode DMComputeL2GradientDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal[], PetscReal *);
447 PETSC_EXTERN PetscErrorCode DMComputeL2FieldDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
448 PETSC_EXTERN PetscErrorCode DMComputeError(DM, Vec, PetscReal[], Vec *);
449 PETSC_EXTERN PetscErrorCode DMHasBasisTransform(DM, PetscBool *);
450 PETSC_EXTERN PetscErrorCode DMCopyTransform(DM, DM);
451 
452 PETSC_EXTERN PetscErrorCode DMGetCompatibility(DM, DM, PetscBool *, PetscBool *);
453 
454 PETSC_EXTERN PetscErrorCode DMMonitorSet(DM, PetscErrorCode (*)(DM, void *), void *, PetscErrorCode (*)(void **));
455 PETSC_EXTERN PetscErrorCode DMMonitorCancel(DM);
456 PETSC_EXTERN PetscErrorCode DMMonitorSetFromOptions(DM, const char[], const char[], const char[], PetscErrorCode (*)(DM, void *), PetscErrorCode (*)(DM, PetscViewerAndFormat *), PetscBool *);
457 PETSC_EXTERN PetscErrorCode DMMonitor(DM);
458 
459 static inline PetscBool DMPolytopeTypeIsHybrid(DMPolytopeType ct)
460 {
461   switch (ct) {
462   case DM_POLYTOPE_POINT_PRISM_TENSOR:
463   case DM_POLYTOPE_SEG_PRISM_TENSOR:
464   case DM_POLYTOPE_TRI_PRISM_TENSOR:
465   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
466     return PETSC_TRUE;
467   default:
468     return PETSC_FALSE;
469   }
470 }
471 
472 static inline PetscInt DMPolytopeTypeGetDim(DMPolytopeType ct)
473 {
474   switch (ct) {
475   case DM_POLYTOPE_POINT:
476     return 0;
477   case DM_POLYTOPE_SEGMENT:
478   case DM_POLYTOPE_POINT_PRISM_TENSOR:
479     return 1;
480   case DM_POLYTOPE_TRIANGLE:
481   case DM_POLYTOPE_QUADRILATERAL:
482   case DM_POLYTOPE_SEG_PRISM_TENSOR:
483   case DM_POLYTOPE_UNKNOWN_FACE:
484     return 2;
485   case DM_POLYTOPE_TETRAHEDRON:
486   case DM_POLYTOPE_HEXAHEDRON:
487   case DM_POLYTOPE_TRI_PRISM:
488   case DM_POLYTOPE_TRI_PRISM_TENSOR:
489   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
490   case DM_POLYTOPE_PYRAMID:
491   case DM_POLYTOPE_UNKNOWN_CELL:
492     return 3;
493   default:
494     return -1;
495   }
496 }
497 
498 static inline PetscInt DMPolytopeTypeGetConeSize(DMPolytopeType ct)
499 {
500   switch (ct) {
501   case DM_POLYTOPE_POINT:
502     return 0;
503   case DM_POLYTOPE_SEGMENT:
504     return 2;
505   case DM_POLYTOPE_POINT_PRISM_TENSOR:
506     return 2;
507   case DM_POLYTOPE_TRIANGLE:
508     return 3;
509   case DM_POLYTOPE_QUADRILATERAL:
510     return 4;
511   case DM_POLYTOPE_SEG_PRISM_TENSOR:
512     return 4;
513   case DM_POLYTOPE_TETRAHEDRON:
514     return 4;
515   case DM_POLYTOPE_HEXAHEDRON:
516     return 6;
517   case DM_POLYTOPE_TRI_PRISM:
518     return 5;
519   case DM_POLYTOPE_TRI_PRISM_TENSOR:
520     return 5;
521   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
522     return 6;
523   case DM_POLYTOPE_PYRAMID:
524     return 5;
525   default:
526     return -1;
527   }
528 }
529 
530 static inline PetscInt DMPolytopeTypeGetNumVertices(DMPolytopeType ct)
531 {
532   switch (ct) {
533   case DM_POLYTOPE_POINT:
534     return 1;
535   case DM_POLYTOPE_SEGMENT:
536     return 2;
537   case DM_POLYTOPE_POINT_PRISM_TENSOR:
538     return 2;
539   case DM_POLYTOPE_TRIANGLE:
540     return 3;
541   case DM_POLYTOPE_QUADRILATERAL:
542     return 4;
543   case DM_POLYTOPE_SEG_PRISM_TENSOR:
544     return 4;
545   case DM_POLYTOPE_TETRAHEDRON:
546     return 4;
547   case DM_POLYTOPE_HEXAHEDRON:
548     return 8;
549   case DM_POLYTOPE_TRI_PRISM:
550     return 6;
551   case DM_POLYTOPE_TRI_PRISM_TENSOR:
552     return 6;
553   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
554     return 8;
555   case DM_POLYTOPE_PYRAMID:
556     return 5;
557   default:
558     return -1;
559   }
560 }
561 
562 static inline DMPolytopeType DMPolytopeTypeSimpleShape(PetscInt dim, PetscBool simplex)
563 {
564   return dim == 0 ? DM_POLYTOPE_POINT : (dim == 1 ? DM_POLYTOPE_SEGMENT : (dim == 2 ? (simplex ? DM_POLYTOPE_TRIANGLE : DM_POLYTOPE_QUADRILATERAL) : (dim == 3 ? (simplex ? DM_POLYTOPE_TETRAHEDRON : DM_POLYTOPE_HEXAHEDRON) : DM_POLYTOPE_UNKNOWN)));
565 }
566 
567 static inline PetscInt DMPolytopeTypeGetNumArrangements(DMPolytopeType ct)
568 {
569   switch (ct) {
570   case DM_POLYTOPE_POINT:
571     return 1;
572   case DM_POLYTOPE_SEGMENT:
573     return 2;
574   case DM_POLYTOPE_POINT_PRISM_TENSOR:
575     return 2;
576   case DM_POLYTOPE_TRIANGLE:
577     return 6;
578   case DM_POLYTOPE_QUADRILATERAL:
579     return 8;
580   case DM_POLYTOPE_SEG_PRISM_TENSOR:
581     return 4;
582   case DM_POLYTOPE_TETRAHEDRON:
583     return 24;
584   case DM_POLYTOPE_HEXAHEDRON:
585     return 48;
586   case DM_POLYTOPE_TRI_PRISM:
587     return 12;
588   case DM_POLYTOPE_TRI_PRISM_TENSOR:
589     return 12;
590   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
591     return 16;
592   case DM_POLYTOPE_PYRAMID:
593     return 8;
594   default:
595     return -1;
596   }
597 }
598 
599 /* An arrangement is a face order combined with an orientation for each face */
600 static inline const PetscInt *DMPolytopeTypeGetArrangement(DMPolytopeType ct, PetscInt o)
601 {
602   static const PetscInt pntArr[1 * 2] = {0, 0};
603   /* a: swap */
604   static const PetscInt segArr[2 * 2 * 2] = {1, 0, 0, 0, /* -1: a */
605                                              0, 0, 1, 0,
606                                              /*  0: e */};
607   /* a: swap first two
608      b: swap last two */
609   static const PetscInt triArr[6 * 3 * 2] = {0, -1, 2, -1, 1, -1, /* -3: b */
610                                              2, -1, 1, -1, 0, -1, /* -2: aba */
611                                              1, -1, 0, -1, 2, -1, /* -1: a */
612                                              0, 0,  1, 0,  2, 0,  /*  0: identity */
613                                              1, 0,  2, 0,  0, 0,  /*  1: ba */
614                                              2, 0,  0, 0,  1, 0,
615                                              /*  2: ab */};
616   /* a: forward cyclic permutation
617      b: swap first and last pairs */
618   static const PetscInt quadArr[8 * 4 * 2] = {1, -1, 0, -1, 3, -1, 2, -1, /* -4: b */
619                                               0, -1, 3, -1, 2, -1, 1, -1, /* -3: b a^3 = a b */
620                                               3, -1, 2, -1, 1, -1, 0, -1, /* -2: b a^2 = a^2 b */
621                                               2, -1, 1, -1, 0, -1, 3, -1, /* -1: b a   = a^3 b */
622                                               0, 0,  1, 0,  2, 0,  3, 0,  /*  0: identity */
623                                               1, 0,  2, 0,  3, 0,  0, 0,  /*  1: a */
624                                               2, 0,  3, 0,  0, 0,  1, 0,  /*  2: a^2 */
625                                               3, 0,  0, 0,  1, 0,  2, 0,
626                                               /*  3: a^3 */};
627   /* r: rotate 180
628      b: swap top and bottom segments */
629   static const PetscInt tsegArr[4 * 4 * 2] = {1, -1, 0, -1, 3, -1, 2, -1, /* -2: r b */
630                                               0, -1, 1, -1, 3, 0,  2, 0,  /* -1: r */
631                                               0, 0,  1, 0,  2, 0,  3, 0,  /*  0: identity */
632                                               1, 0,  0, 0,  2, -1, 3, -1,
633                                               /*  1: b */};
634   /* https://en.wikiversity.org/wiki/Symmetric_group_S4 */
635   static const PetscInt tetArr[24 * 4 * 2] = {3, -2, 2, -3, 0, -1, 1, -1, /* -12: (1324)   p22 */
636                                               3, -1, 1, -3, 2, -1, 0, -1, /* -11: (14)     p21 */
637                                               3, -3, 0, -3, 1, -1, 2, -1, /* -10: (1234)   p18 */
638                                               2, -1, 3, -1, 1, -3, 0, -2, /*  -9: (1423)   p17 */
639                                               2, -3, 0, -1, 3, -2, 1, -3, /*  -8: (1342)   p13 */
640                                               2, -2, 1, -2, 0, -2, 3, -2, /*  -7: (24)     p14 */
641                                               1, -2, 0, -2, 2, -2, 3, -1, /*  -6: (34)     p6  */
642                                               1, -1, 3, -3, 0, -3, 2, -2, /*  -5: (1243)   p10 */
643                                               1, -3, 2, -1, 3, -1, 0, -3, /*  -4: (1432)   p9  */
644                                               0, -3, 1, -1, 3, -3, 2, -3, /*  -3: (12)     p1  */
645                                               0, -2, 2, -2, 1, -2, 3, -3, /*  -2: (23)     p2  */
646                                               0, -1, 3, -2, 2, -3, 1, -2, /*  -1: (13)     p5  */
647                                               0, 0,  1, 0,  2, 0,  3, 0,  /*   0: ()       p0  */
648                                               0, 1,  3, 1,  1, 2,  2, 0,  /*   1: (123)    p4  */
649                                               0, 2,  2, 1,  3, 0,  1, 2,  /*   2: (132)    p3  */
650                                               1, 2,  0, 1,  3, 1,  2, 2,  /*   3: (12)(34) p7  */
651                                               1, 0,  2, 0,  0, 0,  3, 1,  /*   4: (243)    p8  */
652                                               1, 1,  3, 2,  2, 2,  0, 0,  /*   5: (143)    p11 */
653                                               2, 1,  3, 0,  0, 2,  1, 0,  /*   6: (13)(24) p16 */
654                                               2, 2,  1, 1,  3, 2,  0, 2,  /*   7: (142)    p15 */
655                                               2, 0,  0, 0,  1, 0,  3, 2,  /*   8: (234)    p12 */
656                                               3, 2,  2, 2,  1, 1,  0, 1,  /*   9: (14)(23) p23 */
657                                               3, 0,  0, 2,  2, 1,  1, 1,  /*  10: (134)    p19 */
658                                               3, 1,  1, 2,  0, 1,  2, 1 /*  11: (124)    p20 */};
659   /* Each rotation determines a permutation of the four diagonals, and this defines the isomorphism with S_4 */
660   static const PetscInt hexArr[48 * 6 * 2] = {
661     2, -3, 3, -2, 4, -2, 5, -3, 1, -3, 0, -1, /* -24: reflect bottom and use -3 on top */
662     4, -2, 5, -2, 0, -1, 1, -4, 3, -2, 2, -3, /* -23: reflect bottom and use -3 on top */
663     5, -3, 4, -1, 1, -2, 0, -3, 3, -4, 2, -1, /* -22: reflect bottom and use -3 on top */
664     3, -1, 2, -4, 4, -4, 5, -1, 0, -4, 1, -4, /* -21: reflect bottom and use -3 on top */
665     3, -3, 2, -2, 5, -1, 4, -4, 1, -1, 0, -3, /* -20: reflect bottom and use -3 on top */
666     4, -4, 5, -4, 1, -4, 0, -1, 2, -4, 3, -1, /* -19: reflect bottom and use -3 on top */
667     2, -1, 3, -4, 5, -3, 4, -2, 0, -2, 1, -2, /* -18: reflect bottom and use -3 on top */
668     5, -1, 4, -3, 0, -3, 1, -2, 2, -2, 3, -3, /* -17: reflect bottom and use -3 on top */
669     4, -3, 5, -1, 3, -2, 2, -4, 1, -4, 0, -4, /* -16: reflect bottom and use -3 on top */
670     5, -4, 4, -4, 3, -4, 2, -2, 0, -3, 1, -1, /* -15: reflect bottom and use -3 on top */
671     3, -4, 2, -1, 1, -1, 0, -4, 4, -4, 5, -4, /* -14: reflect bottom and use -3 on top */
672     2, -2, 3, -3, 0, -2, 1, -3, 4, -2, 5, -2, /* -13: reflect bottom and use -3 on top */
673     1, -3, 0, -1, 4, -1, 5, -4, 3, -1, 2, -4, /* -12: reflect bottom and use -3 on top */
674     1, -1, 0, -3, 5, -4, 4, -1, 2, -1, 3, -4, /* -11: reflect bottom and use -3 on top */
675     5, -2, 4, -2, 2, -2, 3, -4, 1, -2, 0, -2, /* -10: reflect bottom and use -3 on top */
676     1, -2, 0, -2, 2, -1, 3, -1, 4, -1, 5, -3, /*  -9: reflect bottom and use -3 on top */
677     4, -1, 5, -3, 2, -4, 3, -2, 0, -1, 1, -3, /*  -8: reflect bottom and use -3 on top */
678     3, -2, 2, -3, 0, -4, 1, -1, 5, -1, 4, -3, /*  -7: reflect bottom and use -3 on top */
679     1, -4, 0, -4, 3, -1, 2, -1, 5, -4, 4, -4, /*  -6: reflect bottom and use -3 on top */
680     2, -4, 3, -1, 1, -3, 0, -2, 5, -3, 4, -1, /*  -5: reflect bottom and use -3 on top */
681     0, -4, 1, -4, 4, -3, 5, -2, 2, -3, 3, -2, /*  -4: reflect bottom and use -3 on top */
682     0, -3, 1, -1, 3, -3, 2, -3, 4, -3, 5, -1, /*  -3: reflect bottom and use -3 on top */
683     0, -2, 1, -2, 5, -2, 4, -3, 3, -3, 2, -2, /*  -2: reflect bottom and use -3 on top */
684     0, -1, 1, -3, 2, -3, 3, -3, 5, -2, 4, -2, /*  -1: reflect bottom and use -3 on top */
685     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  5, 0,  /*   0: identity */
686     0, 1,  1, 3,  5, 3,  4, 0,  2, 0,  3, 1,  /*   1: 90  rotation about z */
687     0, 2,  1, 2,  3, 0,  2, 0,  5, 3,  4, 1,  /*   2: 180 rotation about z */
688     0, 3,  1, 1,  4, 0,  5, 3,  3, 0,  2, 1,  /*   3: 270 rotation about z */
689     2, 3,  3, 2,  1, 0,  0, 3,  4, 3,  5, 1,  /*   4: 90  rotation about x */
690     1, 3,  0, 1,  3, 2,  2, 2,  4, 2,  5, 2,  /*   5: 180 rotation about x */
691     3, 1,  2, 0,  0, 1,  1, 2,  4, 1,  5, 3,  /*   6: 270 rotation about x */
692     4, 0,  5, 0,  2, 1,  3, 3,  1, 1,  0, 3,  /*   7: 90  rotation about y */
693     1, 1,  0, 3,  2, 2,  3, 2,  5, 1,  4, 3,  /*   8: 180 rotation about y */
694     5, 1,  4, 3,  2, 3,  3, 1,  0, 0,  1, 0,  /*   9: 270 rotation about y */
695     1, 0,  0, 0,  5, 1,  4, 2,  3, 2,  2, 3,  /*  10: 180 rotation about x+y */
696     1, 2,  0, 2,  4, 2,  5, 1,  2, 2,  3, 3,  /*  11: 180 rotation about x-y */
697     2, 1,  3, 0,  0, 3,  1, 0,  5, 0,  4, 0,  /*  12: 180 rotation about y+z */
698     3, 3,  2, 2,  1, 2,  0, 1,  5, 2,  4, 2,  /*  13: 180 rotation about y-z */
699     5, 3,  4, 1,  3, 1,  2, 3,  1, 3,  0, 1,  /*  14: 180 rotation about z+x */
700     4, 2,  5, 2,  3, 3,  2, 1,  0, 2,  1, 2,  /*  15: 180 rotation about z-x */
701     5, 0,  4, 0,  0, 0,  1, 3,  3, 1,  2, 0,  /*  16: 120 rotation about x+y+z (v0v6) */
702     2, 0,  3, 1,  5, 0,  4, 3,  1, 0,  0, 0,  /*  17: 240 rotation about x+y+z (v0v6) */
703     4, 3,  5, 1,  1, 1,  0, 2,  3, 3,  2, 2,  /*  18: 120 rotation about x+y-z (v4v2) */
704     3, 2,  2, 3,  5, 2,  4, 1,  0, 1,  1, 3,  /*  19: 240 rotation about x+y-z (v4v2) */
705     3, 0,  2, 1,  4, 1,  5, 2,  1, 2,  0, 2,  /*  20: 120 rotation about x-y+z (v1v5) */
706     5, 2,  4, 2,  1, 3,  0, 0,  2, 3,  3, 2,  /*  21: 240 rotation about x-y+z (v1v5) */
707     4, 1,  5, 3,  0, 2,  1, 1,  2, 1,  3, 0,  /*  22: 120 rotation about x-y-z (v7v3) */
708     2, 2,  3, 3,  4, 3,  5, 0,  0, 3,  1, 1,  /*  23: 240 rotation about x-y-z (v7v3) */
709   };
710   static const PetscInt tripArr[12 * 5 * 2] = {
711     1, -3, 0, -1, 3, -1, 4, -1, 2, -1, /* -6: reflect bottom and top */
712     1, -1, 0, -3, 4, -1, 2, -1, 3, -1, /* -5: reflect bottom and top */
713     1, -2, 0, -2, 2, -1, 3, -1, 4, -1, /* -4: reflect bottom and top */
714     0, -3, 1, -1, 3, -3, 2, -3, 4, -3, /* -3: reflect bottom and top */
715     0, -2, 1, -2, 4, -3, 3, -3, 2, -3, /* -2: reflect bottom and top */
716     0, -1, 1, -3, 2, -3, 4, -3, 3, -3, /* -1: reflect bottom and top */
717     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  /*  0: identity */
718     0, 1,  1, 2,  4, 0,  2, 0,  3, 0,  /*  1: 120 rotation about z */
719     0, 2,  1, 1,  3, 0,  4, 0,  2, 0,  /*  2: 240 rotation about z */
720     1, 1,  0, 2,  2, 2,  4, 2,  3, 2,  /*  3: 180 rotation about y of 0 */
721     1, 0,  0, 0,  4, 2,  3, 2,  2, 2,  /*  4: 180 rotation about y of 1 */
722     1, 2,  0, 1,  3, 2,  2, 2,  4, 2,  /*  5: 180 rotation about y of 2 */
723   };
724   /* a: rotate 120 about z
725      b: swap top and bottom segments
726      r: reflect */
727   static const PetscInt ttriArr[12 * 5 * 2] = {
728     1, -3, 0, -3, 2, -2, 4, -2, 3, -2, /* -6: r b a^2 */
729     1, -2, 0, -2, 4, -2, 3, -2, 2, -2, /* -5: r b a */
730     1, -1, 0, -1, 3, -2, 2, -2, 4, -2, /* -4: r b */
731     0, -3, 1, -3, 2, -1, 4, -1, 3, -1, /* -3: r a^2 */
732     0, -2, 1, -2, 4, -1, 3, -1, 2, -1, /* -2: r a */
733     0, -1, 1, -1, 3, -1, 2, -1, 4, -1, /* -1: r */
734     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  /*  0: identity */
735     0, 1,  1, 1,  3, 0,  4, 0,  2, 0,  /*  1: a */
736     0, 2,  1, 2,  4, 0,  2, 0,  3, 0,  /*  2: a^2 */
737     1, 0,  0, 0,  2, 1,  3, 1,  4, 1,  /*  3: b */
738     1, 1,  0, 1,  3, 1,  4, 1,  2, 1,  /*  4: b a */
739     1, 2,  0, 2,  4, 1,  2, 1,  3, 1,  /*  5: b a^2 */
740   };
741   /* a: rotate 90 about z
742      b: swap top and bottom segments
743      r: reflect */
744   static const PetscInt tquadArr[16 * 6 * 2] = {
745     1, -4, 0, -4, 3, -2, 2, -2, 5, -2, 4, -2, /* -8: r b a^3 */
746     1, -3, 0, -3, 2, -2, 5, -2, 4, -2, 3, -2, /* -7: r b a^2 */
747     1, -2, 0, -2, 5, -2, 4, -2, 3, -2, 2, -2, /* -6: r b a */
748     1, -1, 0, -1, 4, -2, 3, -2, 2, -2, 5, -2, /* -5: r b */
749     0, -4, 1, -4, 3, -1, 2, -1, 5, -1, 4, -1, /* -4: r a^3 */
750     0, -3, 1, -3, 2, -1, 5, -1, 4, -1, 3, -1, /* -3: r a^2 */
751     0, -2, 1, -2, 5, -1, 4, -1, 3, -1, 2, -1, /* -2: r a */
752     0, -1, 1, -1, 4, -1, 3, -1, 2, -1, 5, -1, /* -1: r */
753     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  5, 0,  /*  0: identity */
754     0, 1,  1, 1,  3, 0,  4, 0,  5, 0,  2, 0,  /*  1: a */
755     0, 2,  1, 2,  4, 0,  5, 0,  2, 0,  3, 0,  /*  2: a^2 */
756     0, 3,  1, 3,  5, 0,  2, 0,  3, 0,  4, 0,  /*  3: a^3 */
757     1, 0,  0, 0,  2, 1,  3, 1,  4, 1,  5, 1,  /*  4: b */
758     1, 1,  0, 1,  3, 1,  4, 1,  5, 1,  2, 1,  /*  5: b a */
759     1, 2,  0, 2,  4, 1,  5, 1,  2, 1,  3, 1,  /*  6: b a^2 */
760     1, 3,  0, 3,  5, 1,  2, 1,  3, 1,  4, 1,  /*  7: b a^3 */
761   };
762   static const PetscInt pyrArr[8 * 5 * 2] = {
763     0, -4, 2, -3, 1, -3, 4, -3, 3, -3, /* -4: Reflect bottom face */
764     0, -3, 3, -3, 2, -3, 1, -3, 4, -3, /* -3: Reflect bottom face */
765     0, -2, 4, -3, 3, -3, 2, -3, 1, -3, /* -2: Reflect bottom face */
766     0, -1, 1, -3, 4, -3, 3, -3, 2, -3, /* -1: Reflect bottom face */
767     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  /*  0: identity */
768     0, 1,  4, 0,  1, 0,  2, 0,  3, 0,  /*  1:  90 rotation about z */
769     0, 2,  3, 0,  4, 0,  1, 0,  2, 0,  /*  2: 180 rotation about z */
770     0, 3,  2, 0,  3, 0,  4, 0,  1, 0,  /*  3: 270 rotation about z */
771   };
772   switch (ct) {
773   case DM_POLYTOPE_POINT:
774     return pntArr;
775   case DM_POLYTOPE_SEGMENT:
776     return &segArr[(o + 1) * 2 * 2];
777   case DM_POLYTOPE_POINT_PRISM_TENSOR:
778     return &segArr[(o + 1) * 2 * 2];
779   case DM_POLYTOPE_TRIANGLE:
780     return &triArr[(o + 3) * 3 * 2];
781   case DM_POLYTOPE_QUADRILATERAL:
782     return &quadArr[(o + 4) * 4 * 2];
783   case DM_POLYTOPE_SEG_PRISM_TENSOR:
784     return &tsegArr[(o + 2) * 4 * 2];
785   case DM_POLYTOPE_TETRAHEDRON:
786     return &tetArr[(o + 12) * 4 * 2];
787   case DM_POLYTOPE_HEXAHEDRON:
788     return &hexArr[(o + 24) * 6 * 2];
789   case DM_POLYTOPE_TRI_PRISM:
790     return &tripArr[(o + 6) * 5 * 2];
791   case DM_POLYTOPE_TRI_PRISM_TENSOR:
792     return &ttriArr[(o + 6) * 5 * 2];
793   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
794     return &tquadArr[(o + 8) * 6 * 2];
795   case DM_POLYTOPE_PYRAMID:
796     return &pyrArr[(o + 4) * 5 * 2];
797   default:
798     return PETSC_NULLPTR;
799   }
800 }
801 
802 /* A vertex arrangement is a vertex order */
803 static inline const PetscInt *DMPolytopeTypeGetVertexArrangement(DMPolytopeType ct, PetscInt o)
804 {
805   static const PetscInt pntVerts[1]      = {0};
806   static const PetscInt segVerts[2 * 2]  = {1, 0, 0, 1};
807   static const PetscInt triVerts[6 * 3]  = {1, 0, 2, 0, 2, 1, 2, 1, 0, 0, 1, 2, 1, 2, 0, 2, 0, 1};
808   static const PetscInt quadVerts[8 * 4] = {2, 1, 0, 3, 1, 0, 3, 2, 0, 3, 2, 1, 3, 2, 1, 0, 0, 1, 2, 3, 1, 2, 3, 0, 2, 3, 0, 1, 3, 0, 1, 2};
809   static const PetscInt tsegVerts[4 * 4] = {3, 2, 1, 0, 1, 0, 3, 2, 0, 1, 2, 3, 2, 3, 0, 1};
810   static const PetscInt tetVerts[24 * 4] = {2, 3, 1, 0, /* -12: (1324)   p22 */
811                                             3, 1, 2, 0, /* -11: (14)     p21 */
812                                             1, 2, 3, 0, /* -10: (1234)   p18 */
813                                             3, 2, 0, 1, /*  -9: (1423)   p17 */
814                                             2, 0, 3, 1, /*  -8: (1342)   p13 */
815                                             0, 3, 2, 1, /*  -7: (24)     p14 */
816                                             0, 1, 3, 2, /*  -6: (34)     p6  */
817                                             1, 3, 0, 2, /*  -5: (1243)   p10 */
818                                             3, 0, 1, 2, /*  -4: (1432    p9  */
819                                             1, 0, 2, 3, /*  -3: (12)     p1  */
820                                             0, 2, 1, 3, /*  -2: (23)     p2  */
821                                             2, 1, 0, 3, /*  -1: (13)     p5  */
822                                             0, 1, 2, 3, /*   0: ()       p0  */
823                                             1, 2, 0, 3, /*   1: (123)    p4  */
824                                             2, 0, 1, 3, /*   2: (132)    p3  */
825                                             1, 0, 3, 2, /*   3: (12)(34) p7  */
826                                             0, 3, 1, 2, /*   4: (243)    p8  */
827                                             3, 1, 0, 2, /*   5: (143)    p11 */
828                                             2, 3, 0, 1, /*   6: (13)(24) p16 */
829                                             3, 0, 2, 1, /*   7: (142)    p15 */
830                                             0, 2, 3, 1, /*   8: (234)    p12 */
831                                             3, 2, 1, 0, /*   9: (14)(23) p23 */
832                                             2, 1, 3, 0, /*  10: (134)    p19 */
833                                             1, 3, 2, 0 /*  11: (124)    p20 */};
834   static const PetscInt hexVerts[48 * 8] = {
835     3, 0, 4, 5, 2, 6, 7, 1, /* -24: reflected 23 */
836     3, 5, 6, 2, 0, 1, 7, 4, /* -23: reflected 22 */
837     4, 0, 1, 7, 5, 6, 2, 3, /* -22: reflected 21 */
838     6, 7, 1, 2, 5, 3, 0, 4, /* -21: reflected 20 */
839     1, 2, 6, 7, 0, 4, 5, 3, /* -20: reflected 19 */
840     6, 2, 3, 5, 7, 4, 0, 1, /* -19: reflected 18 */
841     4, 5, 3, 0, 7, 1, 2, 6, /* -18: reflected 17 */
842     1, 7, 4, 0, 2, 3, 5, 6, /* -17: reflected 16 */
843     2, 3, 5, 6, 1, 7, 4, 0, /* -16: reflected 15 */
844     7, 4, 0, 1, 6, 2, 3, 5, /* -15: reflected 14 */
845     7, 1, 2, 6, 4, 5, 3, 0, /* -14: reflected 13 */
846     0, 4, 5, 3, 1, 2, 6, 7, /* -13: reflected 12 */
847     5, 4, 7, 6, 3, 2, 1, 0, /* -12: reflected 11 */
848     7, 6, 5, 4, 1, 0, 3, 2, /* -11: reflected 10 */
849     0, 1, 7, 4, 3, 5, 6, 2, /* -10: reflected  9 */
850     4, 7, 6, 5, 0, 3, 2, 1, /*  -9: reflected  8 */
851     5, 6, 2, 3, 4, 0, 1, 7, /*  -8: reflected  7 */
852     2, 6, 7, 1, 3, 0, 4, 5, /*  -7: reflected  6 */
853     6, 5, 4, 7, 2, 1, 0, 3, /*  -6: reflected  5 */
854     5, 3, 0, 4, 6, 7, 1, 2, /*  -5: reflected  4 */
855     2, 1, 0, 3, 6, 5, 4, 7, /*  -4: reflected  3 */
856     1, 0, 3, 2, 7, 6, 5, 4, /*  -3: reflected  2 */
857     0, 3, 2, 1, 4, 7, 6, 5, /*  -2: reflected  1 */
858     3, 2, 1, 0, 5, 4, 7, 6, /*  -1: reflected  0 */
859     0, 1, 2, 3, 4, 5, 6, 7, /*   0: identity */
860     1, 2, 3, 0, 7, 4, 5, 6, /*   1: 90  rotation about z */
861     2, 3, 0, 1, 6, 7, 4, 5, /*   2: 180 rotation about z */
862     3, 0, 1, 2, 5, 6, 7, 4, /*   3: 270 rotation about z */
863     4, 0, 3, 5, 7, 6, 2, 1, /*   4: 90  rotation about x */
864     7, 4, 5, 6, 1, 2, 3, 0, /*   5: 180 rotation about x */
865     1, 7, 6, 2, 0, 3, 5, 4, /*   6: 270 rotation about x */
866     3, 2, 6, 5, 0, 4, 7, 1, /*   7: 90  rotation about y */
867     5, 6, 7, 4, 3, 0, 1, 2, /*   8: 180 rotation about y */
868     4, 7, 1, 0, 5, 3, 2, 6, /*   9: 270 rotation about y */
869     4, 5, 6, 7, 0, 1, 2, 3, /*  10: 180 rotation about x+y */
870     6, 7, 4, 5, 2, 3, 0, 1, /*  11: 180 rotation about x-y */
871     3, 5, 4, 0, 2, 1, 7, 6, /*  12: 180 rotation about y+z */
872     6, 2, 1, 7, 5, 4, 0, 3, /*  13: 180 rotation about y-z */
873     1, 0, 4, 7, 2, 6, 5, 3, /*  14: 180 rotation about z+x */
874     6, 5, 3, 2, 7, 1, 0, 4, /*  15: 180 rotation about z-x */
875     0, 4, 7, 1, 3, 2, 6, 5, /*  16: 120 rotation about x+y+z (v0v6) */
876     0, 3, 5, 4, 1, 7, 6, 2, /*  17: 240 rotation about x+y+z (v0v6) */
877     5, 3, 2, 6, 4, 7, 1, 0, /*  18: 120 rotation about x+y-z (v4v2) */
878     7, 6, 2, 1, 4, 0, 3, 5, /*  19: 240 rotation about x+y-z (v4v2) */
879     2, 1, 7, 6, 3, 5, 4, 0, /*  20: 120 rotation about x-y+z (v1v5) */
880     7, 1, 0, 4, 6, 5, 3, 2, /*  21: 240 rotation about x-y+z (v1v5) */
881     2, 6, 5, 3, 1, 0, 4, 7, /*  22: 120 rotation about x-y-z (v7v3) */
882     5, 4, 0, 3, 6, 2, 1, 7, /*  23: 240 rotation about x-y-z (v7v3) */
883   };
884   static const PetscInt tripVerts[12 * 6] = {
885     4, 3, 5, 2, 1, 0, /* -6: reflect bottom and top */
886     5, 4, 3, 1, 0, 2, /* -5: reflect bottom and top */
887     3, 5, 4, 0, 2, 1, /* -4: reflect bottom and top */
888     1, 0, 2, 5, 4, 3, /* -3: reflect bottom and top */
889     0, 2, 1, 3, 5, 4, /* -2: reflect bottom and top */
890     2, 1, 0, 4, 3, 5, /* -1: reflect bottom and top */
891     0, 1, 2, 3, 4, 5, /*  0: identity */
892     1, 2, 0, 5, 3, 4, /*  1: 120 rotation about z */
893     2, 0, 1, 4, 5, 3, /*  2: 240 rotation about z */
894     4, 5, 3, 2, 0, 1, /*  3: 180 rotation about y of 0 */
895     3, 4, 5, 0, 1, 2, /*  4: 180 rotation about y of 1 */
896     5, 3, 4, 1, 2, 0, /*  5: 180 rotation about y of 2 */
897   };
898   static const PetscInt ttriVerts[12 * 6] = {
899     4, 3, 5, 1, 0, 2, /* -6: r b a^2 */
900     3, 5, 4, 0, 2, 1, /* -5: r b a */
901     5, 4, 3, 2, 1, 0, /* -4: r b */
902     1, 0, 2, 4, 3, 5, /* -3: r a^2 */
903     0, 2, 1, 3, 5, 4, /* -2: r a */
904     2, 1, 0, 5, 4, 3, /* -1: r */
905     0, 1, 2, 3, 4, 5, /*  0: identity */
906     1, 2, 0, 4, 5, 3, /*  1: a */
907     2, 0, 1, 5, 3, 4, /*  2: a^2 */
908     3, 4, 5, 0, 1, 2, /*  3: b */
909     4, 5, 3, 1, 2, 0, /*  4: b a */
910     5, 3, 4, 2, 0, 1, /*  5: b a^2 */
911   };
912   /* a: rotate 90 about z
913      b: swap top and bottom segments
914      r: reflect */
915   static const PetscInt tquadVerts[16 * 8] = {
916     6, 5, 4, 7, 2, 1, 0, 3, /* -8: r b a^3 */
917     5, 4, 7, 6, 1, 0, 3, 2, /* -7: r b a^2 */
918     4, 7, 6, 5, 0, 3, 2, 1, /* -6: r b a */
919     7, 6, 5, 4, 3, 2, 1, 0, /* -5: r b */
920     2, 1, 0, 3, 6, 5, 4, 7, /* -4: r a^3 */
921     1, 0, 3, 2, 5, 4, 7, 6, /* -3: r a^2 */
922     0, 3, 2, 1, 4, 7, 6, 5, /* -2: r a */
923     3, 2, 1, 0, 7, 6, 5, 4, /* -1: r */
924     0, 1, 2, 3, 4, 5, 6, 7, /*  0: identity */
925     1, 2, 3, 0, 5, 6, 7, 4, /*  1: a */
926     2, 3, 0, 1, 6, 7, 4, 5, /*  2: a^2 */
927     3, 0, 1, 2, 7, 4, 5, 6, /*  3: a^3 */
928     4, 5, 6, 7, 0, 1, 2, 3, /*  4: b */
929     5, 6, 7, 4, 1, 2, 3, 0, /*  5: b a */
930     6, 7, 4, 5, 2, 3, 0, 1, /*  6: b a^2 */
931     7, 4, 5, 6, 3, 0, 1, 2, /*  7: b a^3 */
932   };
933   static const PetscInt pyrVerts[8 * 5] = {
934     2, 1, 0, 3, 4, /* -4: Reflect bottom face */
935     1, 0, 3, 2, 4, /* -3: Reflect bottom face */
936     0, 3, 2, 1, 4, /* -2: Reflect bottom face */
937     3, 2, 1, 0, 4, /* -1: Reflect bottom face */
938     0, 1, 2, 3, 4, /*  0: identity */
939     1, 2, 3, 0, 4, /*  1:  90 rotation about z */
940     2, 3, 0, 1, 4, /*  2: 180 rotation about z */
941     3, 0, 1, 2, 4, /*  3: 270 rotation about z */
942   };
943   switch (ct) {
944   case DM_POLYTOPE_POINT:
945     return pntVerts;
946   case DM_POLYTOPE_SEGMENT:
947     return &segVerts[(o + 1) * 2];
948   case DM_POLYTOPE_POINT_PRISM_TENSOR:
949     return &segVerts[(o + 1) * 2];
950   case DM_POLYTOPE_TRIANGLE:
951     return &triVerts[(o + 3) * 3];
952   case DM_POLYTOPE_QUADRILATERAL:
953     return &quadVerts[(o + 4) * 4];
954   case DM_POLYTOPE_SEG_PRISM_TENSOR:
955     return &tsegVerts[(o + 2) * 4];
956   case DM_POLYTOPE_TETRAHEDRON:
957     return &tetVerts[(o + 12) * 4];
958   case DM_POLYTOPE_HEXAHEDRON:
959     return &hexVerts[(o + 24) * 8];
960   case DM_POLYTOPE_TRI_PRISM:
961     return &tripVerts[(o + 6) * 6];
962   case DM_POLYTOPE_TRI_PRISM_TENSOR:
963     return &ttriVerts[(o + 6) * 6];
964   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
965     return &tquadVerts[(o + 8) * 8];
966   case DM_POLYTOPE_PYRAMID:
967     return &pyrVerts[(o + 4) * 5];
968   default:
969     return PETSC_NULLPTR;
970   }
971 }
972 
973 /* This is orientation o1 acting on orientation o2 */
974 static inline PetscInt DMPolytopeTypeComposeOrientation(DMPolytopeType ct, PetscInt o1, PetscInt o2)
975 {
976   static const PetscInt segMult[2 * 2]   = {0, -1, -1, 0};
977   static const PetscInt triMult[6 * 6]   = {0, 2, 1, -3, -1, -2, 1, 0, 2, -2, -3, -1, 2, 1, 0, -1, -2, -3, -3, -2, -1, 0, 1, 2, -2, -1, -3, 1, 2, 0, -1, -3, -2, 2, 0, 1};
978   static const PetscInt quadMult[8 * 8]  = {0,  3,  2,  1,  -4, -1, -2, -3, 1,  0,  3,  2,  -3, -4, -1, -2, 2,  1,  0,  3,  -2, -3, -4, -1, 3,  2,  1,  0,  -1, -2, -3, -4,
979                                             -4, -3, -2, -1, 0,  1,  2,  3,  -3, -2, -1, -4, 1,  2,  3,  0,  -2, -1, -4, -3, 2,  3,  0,  1,  -1, -4, -3, -2, 3,  0,  1,  2};
980   static const PetscInt tsegMult[4 * 4]  = {0, 1, -2, -1, 1, 0, -1, -2, -2, -1, 0, 1, -1, -2, 1, 0};
981   static const PetscInt tetMult[24 * 24] = {
982     3,   2,   7,   0,   5,   10,  9,   8,   1,   6,   11,  4,   -12, -7,  -5,  -9,  -10, -2,  -6,  -1,  -11, -3,  -4,  -8,  4,   0,   8,   1,   3,   11,  10,  6,   2,   7,   9,   5,   -11, -9,  -4,  -8,  -12, -1,  -5,  -3,  -10, -2,  -6,  -7,
983     5,   1,   6,   2,   4,   9,   11,  7,   0,   8,   10,  3,   -10, -8,  -6,  -7,  -11, -3,  -4,  -2,  -12, -1,  -5,  -9,  0,   8,   4,   3,   11,  1,   6,   2,   10,  9,   5,   7,   -9,  -4,  -11, -12, -1,  -8,  -3,  -10, -5,  -6,  -7,  -2,
984     1,   6,   5,   4,   9,   2,   7,   0,   11,  10,  3,   8,   -8,  -6,  -10, -11, -3,  -7,  -2,  -12, -4,  -5,  -9,  -1,  2,   7,   3,   5,   10,  0,   8,   1,   9,   11,  4,   6,   -7,  -5,  -12, -10, -2,  -9,  -1,  -11, -6,  -4,  -8,  -3,
985     6,   5,   1,   9,   2,   4,   0,   11,  7,   3,   8,   10,  -6,  -10, -8,  -3,  -7,  -11, -12, -4,  -2,  -9,  -1,  -5,  7,   3,   2,   10,  0,   5,   1,   9,   8,   4,   6,   11,  -5,  -12, -7,  -2,  -9,  -10, -11, -6,  -1,  -8,  -3,  -4,
986     8,   4,   0,   11,  1,   3,   2,   10,  6,   5,   7,   9,   -4,  -11, -9,  -1,  -8,  -12, -10, -5,  -3,  -7,  -2,  -6,  9,   11,  10,  6,   8,   7,   3,   5,   4,   0,   2,   1,   -3,  -1,  -2,  -6,  -4,  -5,  -9,  -7,  -8,  -12, -10, -11,
987     10,  9,   11,  7,   6,   8,   4,   3,   5,   1,   0,   2,   -2,  -3,  -1,  -5,  -6,  -4,  -8,  -9,  -7,  -11, -12, -10, 11,  10,  9,   8,   7,   6,   5,   4,   3,   2,   1,   0,   -1,  -2,  -3,  -4,  -5,  -6,  -7,  -8,  -9,  -10, -11, -12,
988     -12, -11, -10, -9,  -8,  -7,  -6,  -5,  -4,  -3,  -2,  -1,  0,   1,   2,   3,   4,   5,   6,   7,   8,   9,   10,  11,  -11, -10, -12, -8,  -7,  -9,  -5,  -4,  -6,  -2,  -1,  -3,  1,   2,   0,   4,   5,   3,   7,   8,   6,   10,  11,  9,
989     -10, -12, -11, -7,  -9,  -8,  -4,  -6,  -5,  -1,  -3,  -2,  2,   0,   1,   5,   3,   4,   8,   6,   7,   11,  9,   10,  -9,  -5,  -1,  -12, -2,  -4,  -3,  -11, -7,  -6,  -8,  -10, 3,   10,  8,   0,   7,   11,  9,   4,   2,   6,   1,   5,
990     -8,  -4,  -3,  -11, -1,  -6,  -2,  -10, -9,  -5,  -7,  -12, 4,   11,  6,   1,   8,   9,   10,  5,   0,   7,   2,   3,   -7,  -6,  -2,  -10, -3,  -5,  -1,  -12, -8,  -4,  -9,  -11, 5,   9,   7,   2,   6,   10,  11,  3,   1,   8,   0,   4,
991     -3,  -8,  -4,  -6,  -11, -1,  -9,  -2,  -10, -12, -5,  -7,  6,   4,   11,  9,   1,   8,   0,   10,  5,   3,   7,   2,   -2,  -7,  -6,  -5,  -10, -3,  -8,  -1,  -12, -11, -4,  -9,  7,   5,   9,   10,  2,   6,   1,   11,  3,   4,   8,   0,
992     -1,  -9,  -5,  -4,  -12, -2,  -7,  -3,  -11, -10, -6,  -8,  8,   3,   10,  11,  0,   7,   2,   9,   4,   5,   6,   1,   -6,  -2,  -7,  -3,  -5,  -10, -12, -8,  -1,  -9,  -11, -4,  9,   7,   5,   6,   10,  2,   3,   1,   11,  0,   4,   8,
993     -5,  -1,  -9,  -2,  -4,  -12, -11, -7,  -3,  -8,  -10, -6,  10,  8,   3,   7,   11,  0,   4,   2,   9,   1,   5,   6,   -4,  -3,  -8,  -1,  -6,  -11, -10, -9,  -2,  -7,  -12, -5,  11,  6,   4,   8,   9,   1,   5,   0,   10,  2,   3,   7,
994   };
995   static const PetscInt hexMult[48 * 48] = {
996     18,  2,   5,   22,  21,  8,   16,  0,   13,  6,   11,  3,   15,  9,   4,   23,  12,  1,   19,  10,  7,   20,  14,  17,  -24, -10, -20, -16, -12, -21, -4,  -5,  -18, -13, -15, -8,  -2,  -11, -14, -7,  -3,  -22, -6,  -17, -19, -9,  -1,  -23,
997     8,   20,  19,  2,   5,   23,  0,   17,  11,  1,   15,  7,   13,  4,   10,  18,  3,   14,  21,  9,   12,  22,  6,   16,  -23, -13, -17, -7,  -8,  -19, -16, -12, -22, -2,  -14, -5,  -10, -15, -11, -4,  -20, -9,  -21, -3,  -6,  -18, -24, -1,
998     2,   17,  23,  8,   0,   19,  5,   20,  1,   11,  9,   14,  12,  6,   3,   16,  10,  7,   22,  15,  13,  21,  4,   18,  -22, -14, -19, -5,  -15, -17, -10, -2,  -23, -12, -13, -7,  -16, -8,  -4,  -11, -24, -3,  -18, -9,  -1,  -21, -20, -6,
999     21,  5,   2,   16,  18,  0,   22,  8,   4,   12,  3,   11,  14,  7,   13,  20,  6,   10,  17,  1,   9,   23,  15,  19,  -21, -8,  -18, -15, -4,  -24, -12, -14, -20, -7,  -16, -10, -11, -2,  -5,  -13, -6,  -19, -3,  -23, -22, -1,  -9,  -17,
1000     16,  8,   0,   21,  22,  2,   18,  5,   12,  4,   1,   10,  9,   15,  6,   19,  13,  11,  23,  3,   14,  17,  7,   20,  -20, -16, -24, -10, -2,  -18, -11, -7,  -21, -14, -8,  -15, -12, -4,  -13, -5,  -9,  -23, -1,  -19, -17, -3,  -6,  -22,
1001     5,   19,  20,  0,   8,   17,  2,   23,  10,  3,   7,   15,  6,   12,  11,  22,  1,   9,   16,  14,  4,   18,  13,  21,  -19, -5,  -22, -14, -16, -23, -8,  -11, -17, -4,  -7,  -13, -15, -10, -12, -2,  -21, -6,  -20, -1,  -9,  -24, -18, -3,
1002     22,  0,   8,   18,  16,  5,   21,  2,   6,   13,  10,  1,   7,   14,  12,  17,  4,   3,   20,  11,  15,  19,  9,   23,  -18, -15, -21, -8,  -11, -20, -2,  -13, -24, -5,  -10, -16, -4,  -12, -7,  -14, -1,  -17, -9,  -22, -23, -6,  -3,  -19,
1003     0,   23,  17,  5,   2,   20,  8,   19,  3,   10,  14,  9,   4,   13,  1,   21,  11,  15,  18,  7,   6,   16,  12,  22,  -17, -7,  -23, -13, -10, -22, -15, -4,  -19, -11, -5,  -14, -8,  -16, -2,  -12, -18, -1,  -24, -6,  -3,  -20, -21, -9,
1004     10,  13,  6,   1,   11,  12,  3,   4,   8,   0,   22,  18,  19,  23,  5,   15,  2,   21,  9,   16,  17,  7,   20,  14,  -16, -24, -10, -20, -23, -8,  -19, -6,  -15, -3,  -21, -18, -22, -17, -9,  -1,  -14, -12, -7,  -4,  -11, -13, -5,  -2,
1005     1,   4,   12,  10,  3,   6,   11,  13,  0,   8,   16,  21,  17,  20,  2,   14,  5,   18,  7,   22,  19,  9,   23,  15,  -15, -21, -8,  -18, -17, -10, -22, -3,  -16, -6,  -24, -20, -19, -23, -1,  -9,  -5,  -4,  -13, -12, -2,  -7,  -14, -11,
1006     14,  10,  3,   9,   7,   1,   15,  11,  17,  23,  0,   5,   16,  22,  20,  6,   19,  8,   12,  2,   21,  4,   18,  13,  -14, -19, -5,  -22, -3,  -13, -9,  -20, -7,  -21, -23, -17, -6,  -1,  -24, -18, -12, -16, -2,  -8,  -10, -4,  -11, -15,
1007     7,   3,   10,  15,  14,  11,  9,   1,   20,  19,  5,   0,   18,  21,  17,  4,   23,  2,   13,  8,   22,  6,   16,  12,  -13, -17, -7,  -23, -9,  -14, -3,  -24, -5,  -18, -22, -19, -1,  -6,  -20, -21, -2,  -10, -12, -15, -16, -11, -4,  -8,
1008     13,  14,  15,  12,  4,   9,   6,   7,   21,  22,  23,  20,  2,   0,   18,  3,   16,  17,  1,   19,  8,   11,  5,   10,  -12, -9,  -11, -6,  -21, -4,  -24, -22, -2,  -23, -3,  -1,  -20, -18, -19, -17, -16, -14, -15, -13, -5,  -8,  -10, -7,
1009     6,   9,   7,   4,   12,  14,  13,  15,  16,  18,  17,  19,  0,   2,   22,  1,   21,  23,  3,   20,  5,   10,  8,   11,  -11, -6,  -12, -9,  -20, -2,  -18, -17, -4,  -19, -1,  -3,  -21, -24, -23, -22, -8,  -7,  -10, -5,  -13, -16, -15, -14,
1010     3,   12,  4,   11,  1,   13,  10,  6,   2,   5,   21,  16,  23,  19,  0,   9,   8,   22,  15,  18,  20,  14,  17,  7,   -10, -20, -16, -24, -22, -15, -17, -1,  -8,  -9,  -18, -21, -23, -19, -3,  -6,  -13, -2,  -5,  -11, -4,  -14, -7,  -12,
1011     20,  16,  18,  23,  17,  21,  19,  22,  14,  15,  4,   6,   3,   1,   7,   0,   9,   12,  2,   13,  11,  5,   10,  8,   -9,  -11, -6,  -12, -14, -3,  -13, -10, -1,  -8,  -2,  -4,  -7,  -5,  -16, -15, -23, -20, -22, -18, -24, -19, -17, -21,
1012     11,  6,   13,  3,   10,  4,   1,   12,  5,   2,   18,  22,  20,  17,  8,   7,   0,   16,  14,  21,  23,  15,  19,  9,   -8,  -18, -15, -21, -19, -16, -23, -9,  -10, -1,  -20, -24, -17, -22, -6,  -3,  -7,  -11, -14, -2,  -12, -5,  -13, -4,
1013     9,   11,  1,   14,  15,  3,   7,   10,  23,  17,  2,   8,   21,  18,  19,  13,  20,  5,   4,   0,   16,  12,  22,  6,   -7,  -23, -13, -17, -1,  -5,  -6,  -21, -14, -20, -19, -22, -9,  -3,  -18, -24, -11, -8,  -4,  -16, -15, -2,  -12, -10,
1014     19,  21,  22,  17,  23,  16,  20,  18,  9,   7,   12,  13,  1,   3,   15,  2,   14,  4,   0,   6,   10,  8,   11,  5,   -6,  -12, -9,  -11, -7,  -1,  -5,  -15, -3,  -16, -4,  -2,  -14, -13, -8,  -10, -19, -21, -17, -24, -18, -23, -22, -20,
1015     15,  1,   11,  7,   9,   10,  14,  3,   19,  20,  8,   2,   22,  16,  23,  12,  17,  0,   6,   5,   18,  13,  21,  4,   -5,  -22, -14, -19, -6,  -7,  -1,  -18, -13, -24, -17, -23, -3,  -9,  -21, -20, -4,  -15, -11, -10, -8,  -12, -2,  -16,
1016     4,   15,  14,  6,   13,  7,   12,  9,   18,  16,  20,  23,  5,   8,   21,  11,  22,  19,  10,  17,  0,   3,   2,   1,   -4,  -1,  -2,  -3,  -24, -12, -21, -19, -11, -17, -6,  -9,  -18, -20, -22, -23, -15, -5,  -16, -7,  -14, -10, -8,  -13,
1017     17,  18,  16,  19,  20,  22,  23,  21,  7,   9,   6,   4,   10,  11,  14,  5,   15,  13,  8,   12,  1,   0,   3,   2,   -3,  -4,  -1,  -2,  -13, -9,  -14, -16, -6,  -15, -12, -11, -5,  -7,  -10, -8,  -22, -24, -23, -21, -20, -17, -19, -18,
1018     12,  7,   9,   13,  6,   15,  4,   14,  22,  21,  19,  17,  8,   5,   16,  10,  18,  20,  11,  23,  2,   1,   0,   3,   -2,  -3,  -4,  -1,  -18, -11, -20, -23, -12, -22, -9,  -6,  -24, -21, -17, -19, -10, -13, -8,  -14, -7,  -15, -16, -5,
1019     23,  22,  21,  20,  19,  18,  17,  16,  15,  14,  13,  12,  11,  10,  9,   8,   7,   6,   5,   4,   3,   2,   1,   0,   -1,  -2,  -3,  -4,  -5,  -6,  -7,  -8,  -9,  -10, -11, -12, -13, -14, -15, -16, -17, -18, -19, -20, -21, -22, -23, -24,
1020     -24, -23, -22, -21, -20, -19, -18, -17, -16, -15, -14, -13, -12, -11, -10, -9,  -8,  -7,  -6,  -5,  -4,  -3,  -2,  -1,  0,   1,   2,   3,   4,   5,   6,   7,   8,   9,   10,  11,  12,  13,  14,  15,  16,  17,  18,  19,  20,  21,  22,  23,
1021     -13, -8,  -10, -14, -7,  -16, -5,  -15, -23, -22, -20, -18, -9,  -6,  -17, -11, -19, -21, -12, -24, -3,  -2,  -1,  -4,  1,   2,   3,   0,   17,  10,  19,  22,  11,  21,  8,   5,   23,  20,  16,  18,  9,   12,  7,   13,  6,   14,  15,  4,
1022     -18, -19, -17, -20, -21, -23, -24, -22, -8,  -10, -7,  -5,  -11, -12, -15, -6,  -16, -14, -9,  -13, -2,  -1,  -4,  -3,  2,   3,   0,   1,   12,  8,   13,  15,  5,   14,  11,  10,  4,   6,   9,   7,   21,  23,  22,  20,  19,  16,  18,  17,
1023     -5,  -16, -15, -7,  -14, -8,  -13, -10, -19, -17, -21, -24, -6,  -9,  -22, -12, -23, -20, -11, -18, -1,  -4,  -3,  -2,  3,   0,   1,   2,   23,  11,  20,  18,  10,  16,  5,   8,   17,  19,  21,  22,  14,  4,   15,  6,   13,  9,   7,   12,
1024     -16, -2,  -12, -8,  -10, -11, -15, -4,  -20, -21, -9,  -3,  -23, -17, -24, -13, -18, -1,  -7,  -6,  -19, -14, -22, -5,  4,   21,  13,  18,  5,   6,   0,   17,  12,  23,  16,  22,  2,   8,   20,  19,  3,   14,  10,  9,   7,   11,  1,   15,
1025     -20, -22, -23, -18, -24, -17, -21, -19, -10, -8,  -13, -14, -2,  -4,  -16, -3,  -15, -5,  -1,  -7,  -11, -9,  -12, -6,  5,   11,  8,   10,  6,   0,   4,   14,  2,   15,  3,   1,   13,  12,  7,   9,   18,  20,  16,  23,  17,  22,  21,  19,
1026     -10, -12, -2,  -15, -16, -4,  -8,  -11, -24, -18, -3,  -9,  -22, -19, -20, -14, -21, -6,  -5,  -1,  -17, -13, -23, -7,  6,   22,  12,  16,  0,   4,   5,   20,  13,  19,  18,  21,  8,   2,   17,  23,  10,  7,   3,   15,  14,  1,   11,  9,
1027     -12, -7,  -14, -4,  -11, -5,  -2,  -13, -6,  -3,  -19, -23, -21, -18, -9,  -8,  -1,  -17, -15, -22, -24, -16, -20, -10, 7,   17,  14,  20,  18,  15,  22,  8,   9,   0,   19,  23,  16,  21,  5,   2,   6,   10,  13,  1,   11,  4,   12,  3,
1028     -21, -17, -19, -24, -18, -22, -20, -23, -15, -16, -5,  -7,  -4,  -2,  -8,  -1,  -10, -13, -3,  -14, -12, -6,  -11, -9,  8,   10,  5,   11,  13,  2,   12,  9,   0,   7,   1,   3,   6,   4,   15,  14,  22,  19,  21,  17,  23,  18,  16,  20,
1029     -4,  -13, -5,  -12, -2,  -14, -11, -7,  -3,  -6,  -22, -17, -24, -20, -1,  -10, -9,  -23, -16, -19, -21, -15, -18, -8,  9,   19,  15,  23,  21,  14,  16,  0,   7,   8,   17,  20,  22,  18,  2,   5,   12,  1,   4,   10,  3,   13,  6,   11,
1030     -7,  -10, -8,  -5,  -13, -15, -14, -16, -17, -19, -18, -20, -1,  -3,  -23, -2,  -22, -24, -4,  -21, -6,  -11, -9,  -12, 10,  5,   11,  8,   19,  1,   17,  16,  3,   18,  0,   2,   20,  23,  22,  21,  7,   6,   9,   4,   12,  15,  14,  13,
1031     -14, -15, -16, -13, -5,  -10, -7,  -8,  -22, -23, -24, -21, -3,  -1,  -19, -4,  -17, -18, -2,  -20, -9,  -12, -6,  -11, 11,  8,   10,  5,   20,  3,   23,  21,  1,   22,  2,   0,   19,  17,  18,  16,  15,  13,  14,  12,  4,   7,   9,   6,
1032     -8,  -4,  -11, -16, -15, -12, -10, -2,  -21, -20, -6,  -1,  -19, -22, -18, -5,  -24, -3,  -14, -9,  -23, -7,  -17, -13, 12,  16,  6,   22,  8,   13,  2,   23,  4,   17,  21,  18,  0,   5,   19,  20,  1,   9,   11,  14,  15,  10,  3,   7,
1033     -15, -11, -4,  -10, -8,  -2,  -16, -12, -18, -24, -1,  -6,  -17, -23, -21, -7,  -20, -9,  -13, -3,  -22, -5,  -19, -14, 13,  18,  4,   21,  2,   12,  8,   19,  6,   20,  22,  16,  5,   0,   23,  17,  11,  15,  1,   7,   9,   3,   10,  14,
1034     -2,  -5,  -13, -11, -4,  -7,  -12, -14, -1,  -9,  -17, -22, -18, -21, -3,  -15, -6,  -19, -8,  -23, -20, -10, -24, -16, 14,  20,  7,   17,  16,  9,   21,  2,   15,  5,   23,  19,  18,  22,  0,   8,   4,   3,   12,  11,  1,   6,   13,  10,
1035     -11, -14, -7,  -2,  -12, -13, -4,  -5,  -9,  -1,  -23, -19, -20, -24, -6,  -16, -3,  -22, -10, -17, -18, -8,  -21, -15, 15,  23,  9,   19,  22,  7,   18,  5,   14,  2,   20,  17,  21,  16,  8,   0,   13,  11,  6,   3,   10,  12,  4,   1,
1036     -1,  -24, -18, -6,  -3,  -21, -9,  -20, -4,  -11, -15, -10, -5,  -14, -2,  -22, -12, -16, -19, -8,  -7,  -17, -13, -23, 16,  6,   22,  12,  9,   21,  14,  3,   18,  10,  4,   13,  7,   15,  1,   11,  17,  0,   23,  5,   2,   19,  20,  8,
1037     -23, -1,  -9,  -19, -17, -6,  -22, -3,  -7,  -14, -11, -2,  -8,  -15, -13, -18, -5,  -4,  -21, -12, -16, -20, -10, -24, 17,  14,  20,  7,   10,  19,  1,   12,  23,  4,   9,   15,  3,   11,  6,   13,  0,   16,  8,   21,  22,  5,   2,   18,
1038     -6,  -20, -21, -1,  -9,  -18, -3,  -24, -11, -4,  -8,  -16, -7,  -13, -12, -23, -2,  -10, -17, -15, -5,  -19, -14, -22, 18,  4,   21,  13,  15,  22,  7,   10,  16,  3,   6,   12,  14,  9,   11,  1,   20,  5,   19,  0,   8,   23,  17,  2,
1039     -17, -9,  -1,  -22, -23, -3,  -19, -6,  -13, -5,  -2,  -11, -10, -16, -7,  -20, -14, -12, -24, -4,  -15, -18, -8,  -21, 19,  15,  23,  9,   1,   17,  10,  6,   20,  13,  7,   14,  11,  3,   12,  4,   8,   22,  0,   18,  16,  2,   5,   21,
1040     -22, -6,  -3,  -17, -19, -1,  -23, -9,  -5,  -13, -4,  -12, -15, -8,  -14, -21, -7,  -11, -18, -2,  -10, -24, -16, -20, 20,  7,   17,  14,  3,   23,  11,  13,  19,  6,   15,  9,   10,  1,   4,   12,  5,   18,  2,   22,  21,  0,   8,   16,
1041     -3,  -18, -24, -9,  -1,  -20, -6,  -21, -2,  -12, -10, -15, -13, -7,  -4,  -17, -11, -8,  -23, -16, -14, -22, -5,  -19, 21,  13,  18,  4,   14,  16,  9,   1,   22,  11,  12,  6,   15,  7,   3,   10,  23,  2,   17,  8,   0,   20,  19,  5,
1042     -9,  -21, -20, -3,  -6,  -24, -1,  -18, -12, -2,  -16, -8,  -14, -5,  -11, -19, -4,  -15, -22, -10, -13, -23, -7,  -17, 22,  12,  16,  6,   7,   18,  15,  11,  21,  1,   13,  4,   9,   14,  10,  3,   19,  8,   20,  2,   5,   17,  23,  0,
1043     -19, -3,  -6,  -23, -22, -9,  -17, -1,  -14, -7,  -12, -4,  -16, -10, -5,  -24, -13, -2,  -20, -11, -8,  -21, -15, -18, 23,  9,   19,  15,  11,  20,  3,   4,   17,  12,  14,  7,   1,   10,  13,  6,   2,   21,  5,   16,  18,  8,   0,   22,
1044   };
1045   static const PetscInt tripMult[12 * 12] = {
1046     1,  0,  2,  3,  5,  4,  -6, -4, -5, -2, -3, -1, 0,  2,  1,  4,  3,  5,  -5, -6, -4, -3, -1, -2, 2,  1,  0,  5,  4,  3,  -4, -5, -6, -1, -2, -3, 4,  3,  5,  0,  2,  1,  -3, -1, -2, -5, -6, -4,
1047     3,  5,  4,  1,  0,  2,  -2, -3, -1, -6, -4, -5, 5,  4,  3,  2,  1,  0,  -1, -2, -3, -4, -5, -6, -6, -5, -4, -3, -2, -1, 0,  1,  2,  3,  4,  5,  -4, -6, -5, -2, -1, -3, 1,  2,  0,  5,  3,  4,
1048     -5, -4, -6, -1, -3, -2, 2,  0,  1,  4,  5,  3,  -3, -2, -1, -6, -5, -4, 3,  4,  5,  0,  1,  2,  -1, -3, -2, -5, -4, -6, 4,  5,  3,  2,  0,  1,  -2, -1, -3, -4, -6, -5, 5,  3,  4,  1,  2,  0,
1049   };
1050   static const PetscInt ttriMult[12 * 12] = {
1051     0,  2,  1,  3,  5,  4,  -6, -4, -5, -3, -1, -2, 1,  0,  2,  4,  3,  5,  -5, -6, -4, -2, -3, -1, 2,  1,  0,  5,  4,  3,  -4, -5, -6, -1, -2, -3, 3,  5,  4,  0,  2,  1,  -3, -1, -2, -6, -4, -5,
1052     4,  3,  5,  1,  0,  2,  -2, -3, -1, -5, -6, -4, 5,  4,  3,  2,  1,  0,  -1, -2, -3, -4, -5, -6, -6, -5, -4, -3, -2, -1, 0,  1,  2,  3,  4,  5,  -5, -4, -6, -2, -1, -3, 1,  2,  0,  4,  5,  3,
1053     -4, -6, -5, -1, -3, -2, 2,  0,  1,  5,  3,  4,  -3, -2, -1, -6, -5, -4, 3,  4,  5,  0,  1,  2,  -2, -1, -3, -5, -4, -6, 4,  5,  3,  1,  2,  0,  -1, -3, -2, -4, -6, -5, 5,  3,  4,  2,  0,  1,
1054   };
1055   static const PetscInt tquadMult[16 * 16] = {
1056     0,  3,  2,  1,  4,  7,  6,  5,  -8, -5, -6, -7, -4, -1, -2, -3, 1,  0,  3,  2,  5,  4,  7,  6,  -7, -8, -5, -6, -3, -4, -1, -2, 2,  1,  0,  3,  6,  5,  4,  7,  -6, -7, -8, -5, -2, -3, -4, -1, 3, 2, 1, 0,
1057     7,  6,  5,  4,  -5, -6, -7, -8, -1, -2, -3, -4, 4,  7,  6,  5,  0,  3,  2,  1,  -4, -1, -2, -3, -8, -5, -6, -7, 5,  4,  7,  6,  1,  0,  3,  2,  -3, -4, -1, -2, -7, -8, -5, -6, 6,  5,  4,  7,  2, 1, 0, 3,
1058     -2, -3, -4, -1, -6, -7, -8, -5, 7,  6,  5,  4,  3,  2,  1,  0,  -1, -2, -3, -4, -5, -6, -7, -8, -8, -7, -6, -5, -4, -3, -2, -1, 0,  1,  2,  3,  4,  5,  6,  7,  -7, -6, -5, -8, -3, -2, -1, -4, 1, 2, 3, 0,
1059     5,  6,  7,  4,  -6, -5, -8, -7, -2, -1, -4, -3, 2,  3,  0,  1,  6,  7,  4,  5,  -5, -8, -7, -6, -1, -4, -3, -2, 3,  0,  1,  2,  7,  4,  5,  6,  -4, -3, -2, -1, -8, -7, -6, -5, 4,  5,  6,  7,  0, 1, 2, 3,
1060     -3, -2, -1, -4, -7, -6, -5, -8, 5,  6,  7,  4,  1,  2,  3,  0,  -2, -1, -4, -3, -6, -5, -8, -7, 6,  7,  4,  5,  2,  3,  0,  1,  -1, -4, -3, -2, -5, -8, -7, -6, 7,  4,  5,  6,  3,  0,  1,  2,
1061   };
1062   static const PetscInt pyrMult[8 * 8] = {
1063     0, 3, 2, 1, -4, -1, -2, -3, 1, 0, 3, 2, -3, -4, -1, -2, 2, 1, 0, 3, -2, -3, -4, -1, 3, 2, 1, 0, -1, -2, -3, -4, -4, -3, -2, -1, 0, 1, 2, 3, -3, -2, -1, -4, 1, 2, 3, 0, -2, -1, -4, -3, 2, 3, 0, 1, -1, -4, -3, -2, 3, 0, 1, 2,
1064   };
1065   switch (ct) {
1066   case DM_POLYTOPE_POINT:
1067     return 0;
1068   case DM_POLYTOPE_SEGMENT:
1069   case DM_POLYTOPE_POINT_PRISM_TENSOR:
1070     return segMult[(o1 + 1) * 2 + o2 + 1];
1071   case DM_POLYTOPE_TRIANGLE:
1072     return triMult[(o1 + 3) * 6 + o2 + 3];
1073   case DM_POLYTOPE_QUADRILATERAL:
1074     return quadMult[(o1 + 4) * 8 + o2 + 4];
1075   case DM_POLYTOPE_SEG_PRISM_TENSOR:
1076     return tsegMult[(o1 + 2) * 4 + o2 + 2];
1077   case DM_POLYTOPE_TETRAHEDRON:
1078     return tetMult[(o1 + 12) * 24 + o2 + 12];
1079   case DM_POLYTOPE_HEXAHEDRON:
1080     return hexMult[(o1 + 24) * 48 + o2 + 24];
1081   case DM_POLYTOPE_TRI_PRISM:
1082     return tripMult[(o1 + 6) * 12 + o2 + 6];
1083   case DM_POLYTOPE_TRI_PRISM_TENSOR:
1084     return ttriMult[(o1 + 6) * 12 + o2 + 6];
1085   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
1086     return tquadMult[(o1 + 8) * 16 + o2 + 8];
1087   case DM_POLYTOPE_PYRAMID:
1088     return pyrMult[(o1 + 4) * 8 + o2 + 4];
1089   default:
1090     return 0;
1091   }
1092 }
1093 
1094 /* This is orientation o1 acting on orientation o2^{-1} */
1095 static inline PetscInt DMPolytopeTypeComposeOrientationInv(DMPolytopeType ct, PetscInt o1, PetscInt o2)
1096 {
1097   static const PetscInt triInv[6]    = {-3, -2, -1, 0, 2, 1};
1098   static const PetscInt quadInv[8]   = {-4, -3, -2, -1, 0, 3, 2, 1};
1099   static const PetscInt tetInv[24]   = {-9, -11, -4, -12, -5, -7, -6, -8, -10, -3, -2, -1, 0, 2, 1, 3, 8, 10, 6, 11, 4, 9, 5, 7};
1100   static const PetscInt hexInv[48]   = {-17, -18, -20, -19, -22, -21, -23, -24, -15, -16, -14, -13, -11, -12, -10, -9, -8, -5, -6, -7, -4, -3, -2, -1, 0, 3, 2, 1, 6, 5, 4, 9, 8, 7, 10, 11, 12, 13, 14, 15, 17, 16, 19, 18, 21, 20, 23, 22};
1101   static const PetscInt tripInv[12]  = {-5, -6, -4, -3, -2, -1, 0, 2, 1, 3, 4, 5};
1102   static const PetscInt ttriInv[12]  = {-6, -5, -4, -3, -2, -1, 0, 2, 1, 3, 5, 4};
1103   static const PetscInt tquadInv[16] = {-8, -7, -6, -5, -4, -3, -2, -1, 0, 3, 2, 1, 4, 7, 6, 5};
1104   static const PetscInt pyrInv[8]    = {-4, -3, -2, -1, 0, 3, 2, 1};
1105   switch (ct) {
1106   case DM_POLYTOPE_POINT:
1107     return 0;
1108   case DM_POLYTOPE_SEGMENT:
1109   case DM_POLYTOPE_POINT_PRISM_TENSOR:
1110     return DMPolytopeTypeComposeOrientation(ct, o1, o2);
1111   case DM_POLYTOPE_TRIANGLE:
1112     return DMPolytopeTypeComposeOrientation(ct, o1, triInv[o2 + 3]);
1113   case DM_POLYTOPE_QUADRILATERAL:
1114     return DMPolytopeTypeComposeOrientation(ct, o1, quadInv[o2 + 4]);
1115   case DM_POLYTOPE_SEG_PRISM_TENSOR:
1116     return DMPolytopeTypeComposeOrientation(ct, o1, o2);
1117   case DM_POLYTOPE_TETRAHEDRON:
1118     return DMPolytopeTypeComposeOrientation(ct, o1, tetInv[o2 + 12]);
1119   case DM_POLYTOPE_HEXAHEDRON:
1120     return DMPolytopeTypeComposeOrientation(ct, o1, hexInv[o2 + 24]);
1121   case DM_POLYTOPE_TRI_PRISM:
1122     return DMPolytopeTypeComposeOrientation(ct, o1, tripInv[o2 + 6]);
1123   case DM_POLYTOPE_TRI_PRISM_TENSOR:
1124     return DMPolytopeTypeComposeOrientation(ct, o1, ttriInv[o2 + 6]);
1125   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
1126     return DMPolytopeTypeComposeOrientation(ct, o1, tquadInv[o2 + 8]);
1127   case DM_POLYTOPE_PYRAMID:
1128     return DMPolytopeTypeComposeOrientation(ct, o1, pyrInv[o2 + 4]);
1129   default:
1130     return 0;
1131   }
1132 }
1133 
1134 PETSC_EXTERN PetscErrorCode DMPolytopeMatchOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *);
1135 PETSC_EXTERN PetscErrorCode DMPolytopeMatchVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *);
1136 PETSC_EXTERN PetscErrorCode DMPolytopeGetOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *);
1137 PETSC_EXTERN PetscErrorCode DMPolytopeGetVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *);
1138 PETSC_EXTERN PetscErrorCode DMPolytopeInCellTest(DMPolytopeType, const PetscReal[], PetscBool *);
1139