xref: /petsc/include/petscdm.h (revision 5962854d720b7b8d98c62edd758f00bbb980e600)
1e1589f56SBarry Smith /*
2e1589f56SBarry Smith       Objects to manage the interactions between the mesh data structures and the algebraic objects
3e1589f56SBarry Smith */
4a4963045SJacob Faibussowitsch #pragma once
52c8e378dSBarry Smith #include <petscmat.h>
61e25c274SJed Brown #include <petscdmtypes.h>
7decb47aaSMatthew G. Knepley #include <petscfetypes.h>
82764a2aaSMatthew G. Knepley #include <petscdstypes.h>
9c58f1c22SToby Isaac #include <petscdmlabel.h>
1007218a29SMatthew G. Knepley #include <petscdt.h>
11e1589f56SBarry Smith 
12ac09b921SBarry Smith /* SUBMANSEC = DM */
13ac09b921SBarry Smith 
14607a6623SBarry Smith PETSC_EXTERN PetscErrorCode DMInitializePackage(void);
15e1589f56SBarry Smith 
16014dd563SJed Brown PETSC_EXTERN PetscClassId DM_CLASSID;
17e1589f56SBarry Smith 
18528c63e0SDave May #define DMLOCATEPOINT_POINT_NOT_FOUND -367
19528c63e0SDave May 
2076bdecfbSBarry Smith /*J
2187497f52SBarry Smith     DMType - String with the name of a PETSc `DM`
22e1589f56SBarry Smith 
23e1589f56SBarry Smith    Level: beginner
24e1589f56SBarry Smith 
251cc06b55SBarry Smith .seealso: [](ch_dmbase), `DMSetType()`, `DMCreate()`, `DM`
2676bdecfbSBarry Smith J*/
2719fd82e9SBarry Smith typedef const char *DMType;
28e1589f56SBarry Smith #define DMDA        "da"
29e1589f56SBarry Smith #define DMCOMPOSITE "composite"
30e1589f56SBarry Smith #define DMSLICED    "sliced"
31fe1899a2SJed Brown #define DMSHELL     "shell"
32ab7f58a0SBarry Smith #define DMPLEX      "plex"
338ac4e037SJed Brown #define DMREDUNDANT "redundant"
343a19ef87SMatthew G Knepley #define DMPATCH     "patch"
351d72bce8STim Tautges #define DMMOAB      "moab"
365f2c45f1SShri Abhyankar #define DMNETWORK   "network"
37ef51cf95SToby Isaac #define DMFOREST    "forest"
38db4d5e8cSToby Isaac #define DMP4EST     "p4est"
39db4d5e8cSToby Isaac #define DMP8EST     "p8est"
402fd35b1fSDave May #define DMSWARM     "swarm"
41d852a638SPatrick Sanan #define DMPRODUCT   "product"
42a3101111SPatrick Sanan #define DMSTAG      "stag"
43e1589f56SBarry Smith 
44bff4a2f0SMatthew G. Knepley PETSC_EXTERN const char *const       DMBoundaryTypes[];
4540967b3bSMatthew G. Knepley PETSC_EXTERN const char *const       DMBoundaryConditionTypes[];
46863027abSJed Brown PETSC_EXTERN const char *const       DMBlockingTypes[];
47140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList       DMList;
48c0517cd5SMatthew G. Knepley PETSC_EXTERN DMGeneratorFunctionList DMGenerateList;
49014dd563SJed Brown PETSC_EXTERN PetscErrorCode          DMCreate(MPI_Comm, DM *);
5038221697SMatthew G. Knepley PETSC_EXTERN PetscErrorCode          DMClone(DM, DM *);
5119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode          DMSetType(DM, DMType);
5219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode          DMGetType(DM, DMType *);
53bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode          DMRegister(const char[], PetscErrorCode (*)(DM));
54014dd563SJed Brown PETSC_EXTERN PetscErrorCode          DMRegisterDestroy(void);
55e1589f56SBarry Smith 
56014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMView(DM, PetscViewer);
57014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLoad(DM, PetscViewer);
58014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDestroy(DM *);
59014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM, Vec *);
60014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM, Vec *);
61014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM, Vec *);
62014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM, Vec *);
63014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM, Vec *);
64014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM, Vec *);
65014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
6650eeb1caSToby Isaac PETSC_EXTERN PetscErrorCode DMClearLocalVectors(DM);
67974ca4ecSStefano Zampini PETSC_EXTERN PetscErrorCode DMClearNamedGlobalVectors(DM);
68974ca4ecSStefano Zampini PETSC_EXTERN PetscErrorCode DMClearNamedLocalVectors(DM);
69e77ac854SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasNamedGlobalVector(DM, const char *, PetscBool *);
70014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM, const char *, Vec *);
71014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM, const char *, Vec *);
72e77ac854SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasNamedLocalVector(DM, const char *, PetscBool *);
732348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM, const char *, Vec *);
742348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM, const char *, Vec *);
75014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM, ISLocalToGlobalMapping *);
76014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM, PetscInt *, char ***, IS **);
77014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM, PetscInt *);
78b412c318SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateColoring(DM, ISColoringType, ISColoring *);
79b412c318SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM, Mat *);
80aa0f6e3cSJed Brown PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateSkip(DM, PetscBool);
81014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM, PetscBool);
82b06ff27eSHong Zhang PETSC_EXTERN PetscErrorCode DMSetMatrixStructureOnly(DM, PetscBool);
83863027abSJed Brown PETSC_EXTERN PetscErrorCode DMSetBlockingType(DM, DMBlockingType);
84863027abSJed Brown PETSC_EXTERN PetscErrorCode DMGetBlockingType(DM, DMBlockingType *);
85014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM, DM, Mat *, Vec *);
863ad4599aSBarry Smith PETSC_EXTERN PetscErrorCode DMCreateRestriction(DM, DM, Mat *);
876dbf9973SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMCreateInjection(DM, DM, Mat *);
88bd041c0cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateMassMatrix(DM, DM, Mat *);
89b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateMassMatrixLumped(DM, Vec *);
9069291d52SBarry Smith PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM, PetscInt, MPI_Datatype, void *);
9169291d52SBarry Smith PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM, PetscInt, MPI_Datatype, void *);
92014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefine(DM, MPI_Comm, DM *);
93014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsen(DM, MPI_Comm, DM *);
94a8fb8f29SToby Isaac PETSC_EXTERN PetscErrorCode DMGetCoarseDM(DM, DM *);
95a8fb8f29SToby Isaac PETSC_EXTERN PetscErrorCode DMSetCoarseDM(DM, DM);
9688bdff64SToby Isaac PETSC_EXTERN PetscErrorCode DMGetFineDM(DM, DM *);
9788bdff64SToby Isaac PETSC_EXTERN PetscErrorCode DMSetFineDM(DM, DM);
98014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM, PetscInt, DM[]);
99014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM, PetscInt, DM[]);
100014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, Vec, Mat, DM, void *), void *);
101dc822a44SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, Vec, Mat, DM, void *), void *);
102014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, DM, void *), void *);
1033d8e3701SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, DM, void *), void *);
104014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestrict(DM, Mat, Vec, Mat, DM);
105014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMInterpolate(DM, Mat, DM);
1061f3379b2SToby Isaac PETSC_EXTERN PetscErrorCode DMInterpolateSolution(DM, DM, Mat, Vec, Vec);
107d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMExtrude(DM, PetscInt, DM *);
108014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
109fe2efc57SMark PETSC_EXTERN PetscErrorCode DMViewFromOptions(DM, PetscObject, const char[]);
110ca266f36SBarry Smith 
111c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerate(DM, const char[], PetscBool, DM *);
1129fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMGenerateRegister(const char[], PetscErrorCode (*)(DM, PetscBool, DM *), PetscErrorCode (*)(DM, PetscReal *, DM *), PetscErrorCode (*)(DM, Vec, DMLabel, DMLabel, DM *), PetscInt);
113c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerateRegisterAll(void);
114c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerateRegisterDestroy(void);
115a1b0c543SToby Isaac PETSC_EXTERN PetscErrorCode DMAdaptLabel(DM, DMLabel, DM *);
1169fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMAdaptMetric(DM, Vec, DMLabel, DMLabel, DM *);
117df0b854cSToby Isaac 
118014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetUp(DM);
119014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM, DM, Mat, Vec *);
120edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMDACreateAggregates()", ) PetscErrorCode DMCreateAggregates(DM, DM, Mat *);
121baf369e7SPeter Brune PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM, PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), void *);
122d4d07f1eSToby Isaac PETSC_EXTERN PetscErrorCode DMLocalToGlobalHookAdd(DM, PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), void *);
12301729b5cSPatrick Sanan PETSC_EXTERN PetscErrorCode DMGlobalToLocal(DM, Vec, InsertMode, Vec);
124014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM, Vec, InsertMode, Vec);
125014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM, Vec, InsertMode, Vec);
12601729b5cSPatrick Sanan PETSC_EXTERN PetscErrorCode DMLocalToGlobal(DM, Vec, InsertMode, Vec);
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM, Vec, InsertMode, Vec);
128014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM, Vec, InsertMode, Vec);
129d78e899eSRichard Tran Mills PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM, Vec, InsertMode, Vec);
130d78e899eSRichard Tran Mills PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM, Vec, InsertMode, Vec);
13119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMConvert(DM, DMType, DM *);
132e1589f56SBarry Smith 
133c73cfb54SMatthew G. Knepley /* Topology support */
134c73cfb54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDimension(DM, PetscInt *);
135c73cfb54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetDimension(DM, PetscInt);
136793f3fe5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDimPoints(DM, PetscInt, PetscInt *, PetscInt *);
1378e4ac7eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetUseNatural(DM, PetscBool *);
1388e4ac7eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetUseNatural(DM, PetscBool);
1396858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNeighbors(DM, PetscInt *, const PetscMPIInt **);
14046e270d4SMatthew G. Knepley 
14146e270d4SMatthew G. Knepley /* Coordinate support */
1426636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM, DM *);
1431cfe2091SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM, DM);
1446858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinateDM(DM, DM *);
1456858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinateDM(DM, DM);
14646e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateDim(DM, PetscInt *);
14746e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateDim(DM, PetscInt);
148e8abe2deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM, PetscSection *);
14946e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM, PetscInt, PetscSection);
1506858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinateSection(DM, PetscSection *);
1516858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinateSection(DM, PetscInt, PetscSection);
1526636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM, Vec *);
1536636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM, Vec);
1546858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinates(DM, Vec *);
1556858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinates(DM, Vec);
15681e9a530SVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalSetUp(DM);
1576858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM, Vec *);
15881e9a530SVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalNoncollective(DM, Vec *);
1592db98f8dSVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalTuple(DM, IS, PetscSection *, Vec *);
1606636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM, Vec);
1616858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM);
1626858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocal(DM, Vec *);
1636858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM, Vec *);
1646858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinatesLocal(DM, Vec);
1656858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateField(DM, DMField *);
1666858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateField(DM, DMField);
1676858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetLocalBoundingBox(DM, PetscReal[], PetscReal[]);
1686858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetBoundingBox(DM, PetscReal[], PetscReal[]);
1696858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectCoordinates(DM, PetscFE);
17062a38674SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocatePoints(DM, Vec, DMPointLocationType, PetscSF *);
1716858538eSMatthew G. Knepley 
1726858538eSMatthew G. Knepley /* Periodicity support */
1734fb89dddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM, const PetscReal *[], const PetscReal *[], const PetscReal *[]);
1744fb89dddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM, const PetscReal[], const PetscReal[], const PetscReal[]);
175e907e85cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate(DM, const PetscScalar[], PetscBool, PetscScalar[]);
1762e17dfb7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocalizeCoordinates(DM);
17736447a5eSToby Isaac PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalized(DM, PetscBool *);
1788f700142SStefano Zampini PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalizedLocal(DM, PetscBool *);
1796636e97aSMatthew G Knepley 
1805dbd56e3SPeter Brune /* block hook interface */
181be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, VecScatter, VecScatter, DM, void *), void *);
182b3a6b972SJed Brown PETSC_EXTERN PetscErrorCode DMSubDomainHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, VecScatter, VecScatter, DM, void *), void *);
183be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM, VecScatter, VecScatter, DM);
1845dbd56e3SPeter Brune 
185014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM, const char[]);
18631697293SDave May PETSC_EXTERN PetscErrorCode DMAppendOptionsPrefix(DM, const char[]);
18731697293SDave May PETSC_EXTERN PetscErrorCode DMGetOptionsPrefix(DM, const char *[]);
18819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetVecType(DM, VecType);
189c0dedaeaSBarry Smith PETSC_EXTERN PetscErrorCode DMGetVecType(DM, VecType *);
19019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetMatType(DM, MatType);
191c0dedaeaSBarry Smith PETSC_EXTERN PetscErrorCode DMGetMatType(DM, MatType *);
1928f1509bcSBarry Smith PETSC_EXTERN PetscErrorCode DMSetISColoringType(DM, ISColoringType);
1938f1509bcSBarry Smith PETSC_EXTERN PetscErrorCode DMGetISColoringType(DM, ISColoringType *);
194014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM, void *);
195014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM, PetscErrorCode (*)(void **));
196014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM, void *);
197014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM, PetscErrorCode (*)(DM, Vec, Vec));
198014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM, PetscBool *);
199b0ae01b7SPeter Brune PETSC_EXTERN PetscErrorCode DMHasColoring(DM, PetscBool *);
2003ad4599aSBarry Smith PETSC_EXTERN PetscErrorCode DMHasCreateRestriction(DM, PetscBool *);
201a7058e45SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMHasCreateInjection(DM, PetscBool *);
202014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM, Vec, Vec);
20393d92d96SBarry Smith 
20437bc7515SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, const PetscInt[], IS *, DM *);
2052adcc780SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSuperDM(DM[], PetscInt, IS **, DM *);
206792b654fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSectionSubDM(DM, PetscInt, const PetscInt[], IS *, DM *);
207792b654fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSectionSuperDM(DM[], PetscInt, IS **, DM *);
20816621825SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM, PetscInt *, char ***, IS **, DM **);
2098d4ac253SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM, PetscInt *, char ***, IS **, IS **, DM **);
210e30e807fSPeter Brune PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM, PetscInt, DM *, VecScatter **, VecScatter **, VecScatter **);
211e7c4fc90SDmitry Karpeev 
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM, PetscInt *);
213fef3a512SBarry Smith PETSC_EXTERN PetscErrorCode DMSetRefineLevel(DM, PetscInt);
214014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM, PetscInt *);
2159a64c4a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoarsenLevel(DM, PetscInt);
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);
217e1589f56SBarry Smith 
2185f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM *);
2195f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
220c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM *);
221c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
222531c7667SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringUseDM(Mat, MatFDColoring);
2235f1ad066SMatthew G Knepley 
224e1589f56SBarry Smith typedef struct NLF_DAAD *NLF;
225e1589f56SBarry Smith 
226bc2bf880SBarry Smith #define DM_FILE_CLASSID 1211221
2277da65231SMatthew G Knepley 
2287da65231SMatthew G Knepley /* FEM support */
2292ecf6ec3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPrintCellIndices(PetscInt, const char[], PetscInt, const PetscInt[]);
230014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char[], PetscInt, const PetscScalar[]);
231*5962854dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPrintCellVectorReal(PetscInt, const char[], PetscInt, const PetscReal[]);
232014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char[], PetscInt, PetscInt, const PetscScalar[]);
2336113b454SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char[], PetscReal, Vec);
2347da65231SMatthew G Knepley 
2358cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *));
2368cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *));
2378cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *));
2388cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *));
239f9d4088aSMatthew G. Knepley 
240061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMGetSection(DM, PetscSection *); /* Use DMGetLocalSection() in new code (since v3.12) */
241061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMSetSection(DM, PetscSection);   /* Use DMSetLocalSection() in new code (since v3.12) */
242061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalSection(DM, PetscSection *);
243061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMSetLocalSection(DM, PetscSection);
244e87a4003SBarry Smith PETSC_EXTERN PetscErrorCode DMGetGlobalSection(DM, PetscSection *);
245e87a4003SBarry Smith PETSC_EXTERN PetscErrorCode DMSetGlobalSection(DM, PetscSection);
246d2b2dc1eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMUseTensorOrder(DM, PetscBool);
247edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMGetSection()", ) PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *s)
248d71ae5a4SJacob Faibussowitsch {
2499371c9d4SSatish Balay   return DMGetSection(dm, s);
2509371c9d4SSatish Balay }
251edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMSetSection()", ) PetscErrorCode DMSetDefaultSection(DM dm, PetscSection s)
252d71ae5a4SJacob Faibussowitsch {
2539371c9d4SSatish Balay   return DMSetSection(dm, s);
2549371c9d4SSatish Balay }
255edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMGetGlobalSection()", ) PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *s)
256d71ae5a4SJacob Faibussowitsch {
2579371c9d4SSatish Balay   return DMGetGlobalSection(dm, s);
2589371c9d4SSatish Balay }
259edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMSetGlobalSection()", ) PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection s)
260d71ae5a4SJacob Faibussowitsch {
2619371c9d4SSatish Balay   return DMSetGlobalSection(dm, s);
2629371c9d4SSatish Balay }
263e87a4003SBarry Smith 
2641bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMGetSectionSF(DM, PetscSF *);
2651bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMSetSectionSF(DM, PetscSF);
2661bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateSectionSF(DM, PetscSection, PetscSection);
267edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMGetSectionSF()", ) PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *s)
268d71ae5a4SJacob Faibussowitsch {
2699371c9d4SSatish Balay   return DMGetSectionSF(dm, s);
2709371c9d4SSatish Balay }
271edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMSetSectionSF()", ) PetscErrorCode DMSetDefaultSF(DM dm, PetscSF s)
272d71ae5a4SJacob Faibussowitsch {
2739371c9d4SSatish Balay   return DMSetSectionSF(dm, s);
2749371c9d4SSatish Balay }
275edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMCreateSectionSF()", ) PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection l, PetscSection g)
276d71ae5a4SJacob Faibussowitsch {
2779371c9d4SSatish Balay   return DMCreateSectionSF(dm, l, g);
2789371c9d4SSatish Balay }
2791bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
2801bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);
2814f37162bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNaturalSF(DM, PetscSF *);
2824f37162bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNaturalSF(DM, PetscSF);
2831bb6d2a8SBarry Smith 
28479769bd5SJed Brown PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *, Vec *);
28579769bd5SJed Brown PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat, Vec);
286e87a4003SBarry Smith 
28714f150ffSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *);
288cdb7a50dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *);
289cdb7a50dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal);
290cdb7a50dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char *, PetscInt, PetscReal *);
29114f150ffSMatthew G. Knepley 
292af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
293af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
29444a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, DMLabel *, PetscObject *);
29544a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, DMLabel, PetscObject);
29644a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAddField(DM, DMLabel, PetscObject);
297e0b68406SMatthew Knepley PETSC_EXTERN PetscErrorCode DMSetFieldAvoidTensor(DM, PetscInt, PetscBool);
298e0b68406SMatthew Knepley PETSC_EXTERN PetscErrorCode DMGetFieldAvoidTensor(DM, PetscInt, PetscBool *);
29944a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMClearFields(DM);
300e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyFields(DM, DM);
30134aa8a36SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAdjacency(DM, PetscInt, PetscBool *, PetscBool *);
30234aa8a36SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetAdjacency(DM, PetscInt, PetscBool, PetscBool);
303b0441da4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetBasicAdjacency(DM, PetscBool *, PetscBool *);
304b0441da4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetBasicAdjacency(DM, PetscBool, PetscBool);
305e5e52638SMatthew G. Knepley 
306e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNumDS(DM, PetscInt *);
307e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *);
30807218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellDS(DM, PetscInt, PetscDS *, PetscDS *);
30907218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetRegionDS(DM, DMLabel, IS *, PetscDS *, PetscDS *);
31007218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetRegionDS(DM, DMLabel, IS, PetscDS, PetscDS);
31107218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetRegionNumDS(DM, PetscInt, DMLabel *, IS *, PetscDS *, PetscDS *);
31207218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetRegionNumDS(DM, PetscInt, DMLabel, IS, PetscDS, PetscDS);
3131d3af9e0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMFindRegionNum(DM, PetscDS, PetscInt *);
3142df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateFEDefault(DM, PetscInt, const char[], PetscInt, PetscFE *);
315e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateDS(DM);
316e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMClearDS(DM);
317e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyDS(DM, DM);
318e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyDisc(DM, DM);
319f2cacb80SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMComputeExactSolution(DM, PetscReal, Vec, Vec);
3209a2a23afSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNumAuxiliaryVec(DM, PetscInt *);
321ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAuxiliaryVec(DM, DMLabel, PetscInt, PetscInt, Vec *);
322ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetAuxiliaryVec(DM, DMLabel, PetscInt, PetscInt, Vec);
323ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAuxiliaryLabels(DM, DMLabel[], PetscInt[], PetscInt[]);
3249a2a23afSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyAuxiliaryVec(DM, DM);
325af122d2aSMatthew G Knepley 
3264267b1a3SMatthew G. Knepley /*MC
3274267b1a3SMatthew G. Knepley   DMInterpolationInfo - Structure for holding information about interpolation on a mesh
3284267b1a3SMatthew G. Knepley 
3294267b1a3SMatthew G. Knepley   Synopsis:
3304267b1a3SMatthew G. Knepley     comm   - The communicator
3314267b1a3SMatthew G. Knepley     dim    - The spatial dimension of points
3324267b1a3SMatthew G. Knepley     nInput - The number of input points
3334267b1a3SMatthew G. Knepley     points - The input point coordinates
3344267b1a3SMatthew G. Knepley     cells  - The cell containing each point
3354267b1a3SMatthew G. Knepley     n      - The number of local points
3364267b1a3SMatthew G. Knepley     coords - The point coordinates
3374267b1a3SMatthew G. Knepley     dof    - The number of components to interpolate
3384267b1a3SMatthew G. Knepley 
33916a05f60SBarry Smith   Level: intermediate
34016a05f60SBarry Smith 
34116a05f60SBarry Smith .seealso: `DM`, `DMInterpolationCreate()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
3424267b1a3SMatthew G. Knepley M*/
343e87bb0d3SMatthew G Knepley struct _DMInterpolationInfo {
344e87bb0d3SMatthew G Knepley   MPI_Comm   comm;
345e87bb0d3SMatthew G Knepley   PetscInt   dim;    /*1 The spatial dimension of points */
346e87bb0d3SMatthew G Knepley   PetscInt   nInput; /* The number of input points */
347e87bb0d3SMatthew G Knepley   PetscReal *points; /* The input point coordinates */
348e87bb0d3SMatthew G Knepley   PetscInt  *cells;  /* The cell containing each point */
349e87bb0d3SMatthew G Knepley   PetscInt   n;      /* The number of local points */
350e87bb0d3SMatthew G Knepley   Vec        coords; /* The point coordinates */
351e87bb0d3SMatthew G Knepley   PetscInt   dof;    /* The number of components to interpolate */
352e87bb0d3SMatthew G Knepley };
353e87bb0d3SMatthew G Knepley typedef struct _DMInterpolationInfo *DMInterpolationInfo;
354e87bb0d3SMatthew G Knepley 
35594b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
35694b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
35794b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
35894b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
35994b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
36094b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
36152aa1562SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool, PetscBool);
36294b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
36394b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
36494b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
36594b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
36694b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);
367c58f1c22SToby Isaac 
368c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMCreateLabel(DM, const char[]);
3694bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreateLabelAtIndex(DM, PetscInt, const char[]);
370c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelValue(DM, const char[], PetscInt, PetscInt *);
371c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMSetLabelValue(DM, const char[], PetscInt, PetscInt);
372c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMClearLabelValue(DM, const char[], PetscInt, PetscInt);
373c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelSize(DM, const char[], PetscInt *);
374c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelIdIS(DM, const char[], IS *);
375c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetStratumSize(DM, const char[], PetscInt, PetscInt *);
376c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetStratumIS(DM, const char[], PetscInt, IS *);
3774de306b1SToby Isaac PETSC_EXTERN PetscErrorCode DMSetStratumIS(DM, const char[], PetscInt, IS);
378c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMClearLabelStratum(DM, const char[], PetscInt);
379c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelOutput(DM, const char[], PetscBool *);
380c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMSetLabelOutput(DM, const char[], PetscBool);
3815f15299fSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMGetFirstLabeledPoint(DM, DM, DMLabel, PetscInt, const PetscInt *, PetscInt, PetscInt *, PetscDS *);
382c58f1c22SToby Isaac 
3832cbb9b06SVaclav Hapla /*E
38487497f52SBarry Smith    DMCopyLabelsMode - Determines how `DMCopyLabels()` behaves when there is a `DMLabel` in the source and destination DMs with the same name
3852cbb9b06SVaclav Hapla 
38616a05f60SBarry Smith    Values:
38716a05f60SBarry Smith +  `DM_COPY_LABELS_REPLACE`  - replace label in destination by label from source
38816a05f60SBarry Smith .  `DM_COPY_LABELS_KEEP`     - keep destination label
38916a05f60SBarry Smith -  `DM_COPY_LABELS_FAIL`     - throw error
39016a05f60SBarry Smith 
3912cbb9b06SVaclav Hapla    Level: advanced
3922cbb9b06SVaclav Hapla 
39316a05f60SBarry Smith .seealso: `DMLabel`, `DM`, `DMCompareLabels()`, `DMRemoveLabel()`
3942cbb9b06SVaclav Hapla E*/
3959371c9d4SSatish Balay typedef enum {
3969371c9d4SSatish Balay   DM_COPY_LABELS_REPLACE,
3979371c9d4SSatish Balay   DM_COPY_LABELS_KEEP,
3989371c9d4SSatish Balay   DM_COPY_LABELS_FAIL
3999371c9d4SSatish Balay } DMCopyLabelsMode;
4002cbb9b06SVaclav Hapla PETSC_EXTERN const char *const DMCopyLabelsModes[];
4012cbb9b06SVaclav Hapla 
402c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetNumLabels(DM, PetscInt *);
403c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelName(DM, PetscInt, const char **);
404c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMHasLabel(DM, const char[], PetscBool *);
405c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabel(DM, const char *, DMLabel *);
4064a7ee7d0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetLabel(DM, DMLabel);
407c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelByNum(DM, PetscInt, DMLabel *);
408c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMAddLabel(DM, DMLabel);
409c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMRemoveLabel(DM, const char[], DMLabel *);
410306894acSVaclav Hapla PETSC_EXTERN PetscErrorCode DMRemoveLabelBySelf(DM, DMLabel *, PetscBool);
4112cbb9b06SVaclav Hapla PETSC_EXTERN PetscErrorCode DMCopyLabels(DM, DM, PetscCopyMode, PetscBool, DMCopyLabelsMode emode);
412609dae6eSVaclav Hapla PETSC_EXTERN PetscErrorCode DMCompareLabels(DM, DM, PetscBool *, char **);
413c58f1c22SToby Isaac 
41445480ffeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAddBoundary(DM, DMBoundaryConditionType, const char[], DMLabel, PetscInt, const PetscInt[], PetscInt, PetscInt, const PetscInt[], void (*)(void), void (*)(void), void *, PetscInt *);
415a6ba4734SToby Isaac PETSC_EXTERN PetscErrorCode DMIsBoundaryPoint(DM, PetscInt, PetscBool *);
4164d6f44ffSToby Isaac 
4170709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMProjectFunction(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
4180709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMProjectFunctionLocal(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
4192c53366bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectFunctionLabel(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
4201c531cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectFunctionLabelLocal(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
421191494d9SMatthew G. Knepley 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);
422d29d7c6eSMatthew G. Knepley 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);
4231c531cf8SMatthew G. Knepley 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);
424ece3a9fcSMatthew G. Knepley 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);
4250709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2Diff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
426b698f381SToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2GradientDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal[], PetscReal *);
4271189c1efSToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2FieldDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
4282e4af2aeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMComputeError(DM, Vec, PetscReal[], Vec *);
429ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasBasisTransform(DM, PetscBool *);
430ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyTransform(DM, DM);
4318320bc6fSPatrick Sanan 
4328320bc6fSPatrick Sanan PETSC_EXTERN PetscErrorCode DMGetCompatibility(DM, DM, PetscBool *, PetscBool *);
433c0f0dcc3SMatthew G. Knepley 
434c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitorSet(DM, PetscErrorCode (*)(DM, void *), void *, PetscErrorCode (*)(void **));
435c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitorCancel(DM);
436c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitorSetFromOptions(DM, const char[], const char[], const char[], PetscErrorCode (*)(DM, void *), PetscErrorCode (*)(DM, PetscViewerAndFormat *), PetscBool *);
437c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitor(DM);
438c0f0dcc3SMatthew G. Knepley 
439d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetDim(DMPolytopeType ct)
440d71ae5a4SJacob Faibussowitsch {
441412e9a14SMatthew G. Knepley   switch (ct) {
442d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
443d71ae5a4SJacob Faibussowitsch     return 0;
444412e9a14SMatthew G. Knepley   case DM_POLYTOPE_SEGMENT:
445d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
446d71ae5a4SJacob Faibussowitsch     return 1;
447412e9a14SMatthew G. Knepley   case DM_POLYTOPE_TRIANGLE:
448412e9a14SMatthew G. Knepley   case DM_POLYTOPE_QUADRILATERAL:
449d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
450d71ae5a4SJacob Faibussowitsch     return 2;
451412e9a14SMatthew G. Knepley   case DM_POLYTOPE_TETRAHEDRON:
452412e9a14SMatthew G. Knepley   case DM_POLYTOPE_HEXAHEDRON:
453412e9a14SMatthew G. Knepley   case DM_POLYTOPE_TRI_PRISM:
454412e9a14SMatthew G. Knepley   case DM_POLYTOPE_TRI_PRISM_TENSOR:
455412e9a14SMatthew G. Knepley   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
456d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
457d71ae5a4SJacob Faibussowitsch     return 3;
458d71ae5a4SJacob Faibussowitsch   default:
459d71ae5a4SJacob Faibussowitsch     return -1;
460412e9a14SMatthew G. Knepley   }
461412e9a14SMatthew G. Knepley }
462412e9a14SMatthew G. Knepley 
463d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetConeSize(DMPolytopeType ct)
464d71ae5a4SJacob Faibussowitsch {
465412e9a14SMatthew G. Knepley   switch (ct) {
466d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
467d71ae5a4SJacob Faibussowitsch     return 0;
468d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
469d71ae5a4SJacob Faibussowitsch     return 2;
470d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
471d71ae5a4SJacob Faibussowitsch     return 2;
472d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
473d71ae5a4SJacob Faibussowitsch     return 3;
474d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
475d71ae5a4SJacob Faibussowitsch     return 4;
476d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
477d71ae5a4SJacob Faibussowitsch     return 4;
478d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
479d71ae5a4SJacob Faibussowitsch     return 4;
480d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
481d71ae5a4SJacob Faibussowitsch     return 6;
482d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
483d71ae5a4SJacob Faibussowitsch     return 5;
484d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
485d71ae5a4SJacob Faibussowitsch     return 5;
486d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
487d71ae5a4SJacob Faibussowitsch     return 6;
488d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
489d71ae5a4SJacob Faibussowitsch     return 5;
490d71ae5a4SJacob Faibussowitsch   default:
491d71ae5a4SJacob Faibussowitsch     return -1;
492412e9a14SMatthew G. Knepley   }
493412e9a14SMatthew G. Knepley }
494412e9a14SMatthew G. Knepley 
495d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetNumVertices(DMPolytopeType ct)
496d71ae5a4SJacob Faibussowitsch {
497412e9a14SMatthew G. Knepley   switch (ct) {
498d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
499d71ae5a4SJacob Faibussowitsch     return 1;
500d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
501d71ae5a4SJacob Faibussowitsch     return 2;
502d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
503d71ae5a4SJacob Faibussowitsch     return 2;
504d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
505d71ae5a4SJacob Faibussowitsch     return 3;
506d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
507d71ae5a4SJacob Faibussowitsch     return 4;
508d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
509d71ae5a4SJacob Faibussowitsch     return 4;
510d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
511d71ae5a4SJacob Faibussowitsch     return 4;
512d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
513d71ae5a4SJacob Faibussowitsch     return 8;
514d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
515d71ae5a4SJacob Faibussowitsch     return 6;
516d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
517d71ae5a4SJacob Faibussowitsch     return 6;
518d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
519d71ae5a4SJacob Faibussowitsch     return 8;
520d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
521d71ae5a4SJacob Faibussowitsch     return 5;
522d71ae5a4SJacob Faibussowitsch   default:
523d71ae5a4SJacob Faibussowitsch     return -1;
524412e9a14SMatthew G. Knepley   }
525412e9a14SMatthew G. Knepley }
526412e9a14SMatthew G. Knepley 
527d71ae5a4SJacob Faibussowitsch static inline DMPolytopeType DMPolytopeTypeSimpleShape(PetscInt dim, PetscBool simplex)
528d71ae5a4SJacob Faibussowitsch {
5299371c9d4SSatish Balay   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)));
5309318fe57SMatthew G. Knepley }
5319318fe57SMatthew G. Knepley 
532d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetNumArrangments(DMPolytopeType ct)
533d71ae5a4SJacob Faibussowitsch {
534b5a892a1SMatthew G. Knepley   switch (ct) {
535d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
536d71ae5a4SJacob Faibussowitsch     return 1;
537d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
538d71ae5a4SJacob Faibussowitsch     return 2;
539d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
540d71ae5a4SJacob Faibussowitsch     return 2;
541d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
542d71ae5a4SJacob Faibussowitsch     return 6;
543d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
544d71ae5a4SJacob Faibussowitsch     return 8;
545d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
546d71ae5a4SJacob Faibussowitsch     return 4;
547d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
548d71ae5a4SJacob Faibussowitsch     return 24;
549d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
550d71ae5a4SJacob Faibussowitsch     return 48;
551d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
552d71ae5a4SJacob Faibussowitsch     return 12;
553d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
554d71ae5a4SJacob Faibussowitsch     return 12;
555d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
556d71ae5a4SJacob Faibussowitsch     return 16;
557d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
558d71ae5a4SJacob Faibussowitsch     return 8;
559d71ae5a4SJacob Faibussowitsch   default:
560d71ae5a4SJacob Faibussowitsch     return -1;
561b5a892a1SMatthew G. Knepley   }
562b5a892a1SMatthew G. Knepley }
563b5a892a1SMatthew G. Knepley 
564b5a892a1SMatthew G. Knepley /* An arrangement is a face order combined with an orientation for each face */
565d71ae5a4SJacob Faibussowitsch static inline const PetscInt *DMPolytopeTypeGetArrangment(DMPolytopeType ct, PetscInt o)
566d71ae5a4SJacob Faibussowitsch {
567ef8b56bfSJed Brown   static const PetscInt pntArr[1 * 2] = {0, 0};
568b5a892a1SMatthew G. Knepley   /* a: swap */
5699371c9d4SSatish Balay   static const PetscInt segArr[2 * 2 * 2] = {1, 0, 0, 0, /* -1: a */
5709371c9d4SSatish Balay                                              0, 0, 1, 0,
5719371c9d4SSatish Balay                                              /*  0: e */};
572b5a892a1SMatthew G. Knepley   /* a: swap first two
573b5a892a1SMatthew G. Knepley      b: swap last two */
5749371c9d4SSatish Balay   static const PetscInt triArr[6 * 3 * 2] = {0, -1, 2, -1, 1, -1, /* -3: b */
575b5a892a1SMatthew G. Knepley                                              2, -1, 1, -1, 0, -1, /* -2: aba */
576b5a892a1SMatthew G. Knepley                                              1, -1, 0, -1, 2, -1, /* -1: a */
577b5a892a1SMatthew G. Knepley                                              0, 0,  1, 0,  2, 0,  /*  0: identity */
578b5a892a1SMatthew G. Knepley                                              1, 0,  2, 0,  0, 0,  /*  1: ba */
5799371c9d4SSatish Balay                                              2, 0,  0, 0,  1, 0,
5809371c9d4SSatish Balay                                              /*  2: ab */};
581b5a892a1SMatthew G. Knepley   /* a: forward cyclic permutation
582b5a892a1SMatthew G. Knepley      b: swap first and last pairs */
5839371c9d4SSatish Balay   static const PetscInt quadArr[8 * 4 * 2] = {1, -1, 0, -1, 3, -1, 2, -1, /* -4: b */
584b5a892a1SMatthew G. Knepley                                               0, -1, 3, -1, 2, -1, 1, -1, /* -3: b a^3 = a b */
585b5a892a1SMatthew G. Knepley                                               3, -1, 2, -1, 1, -1, 0, -1, /* -2: b a^2 = a^2 b */
586b5a892a1SMatthew G. Knepley                                               2, -1, 1, -1, 0, -1, 3, -1, /* -1: b a   = a^3 b */
587b5a892a1SMatthew G. Knepley                                               0, 0,  1, 0,  2, 0,  3, 0,  /*  0: identity */
588b5a892a1SMatthew G. Knepley                                               1, 0,  2, 0,  3, 0,  0, 0,  /*  1: a */
589b5a892a1SMatthew G. Knepley                                               2, 0,  3, 0,  0, 0,  1, 0,  /*  2: a^2 */
5909371c9d4SSatish Balay                                               3, 0,  0, 0,  1, 0,  2, 0,
5919371c9d4SSatish Balay                                               /*  3: a^3 */};
592b5a892a1SMatthew G. Knepley   /* r: rotate 180
593b5a892a1SMatthew G. Knepley      b: swap top and bottom segments */
5949371c9d4SSatish Balay   static const PetscInt tsegArr[4 * 4 * 2] = {1, -1, 0, -1, 3, -1, 2, -1, /* -2: r b */
595b5a892a1SMatthew G. Knepley                                               0, -1, 1, -1, 3, 0,  2, 0,  /* -1: r */
596b5a892a1SMatthew G. Knepley                                               0, 0,  1, 0,  2, 0,  3, 0,  /*  0: identity */
5979371c9d4SSatish Balay                                               1, 0,  0, 0,  2, -1, 3, -1,
5989371c9d4SSatish Balay                                               /*  1: b */};
599b5a892a1SMatthew G. Knepley   /* https://en.wikiversity.org/wiki/Symmetric_group_S4 */
6009371c9d4SSatish Balay   static const PetscInt tetArr[24 * 4 * 2] = {3, -2, 2, -3, 0, -1, 1, -1, /* -12: (1324)   p22 */
601b5a892a1SMatthew G. Knepley                                               3, -1, 1, -3, 2, -1, 0, -1, /* -11: (14)     p21 */
602b5a892a1SMatthew G. Knepley                                               3, -3, 0, -3, 1, -1, 2, -1, /* -10: (1234)   p18 */
603b5a892a1SMatthew G. Knepley                                               2, -1, 3, -1, 1, -3, 0, -2, /*  -9: (1423)   p17 */
604b5a892a1SMatthew G. Knepley                                               2, -3, 0, -1, 3, -2, 1, -3, /*  -8: (1342)   p13 */
605b5a892a1SMatthew G. Knepley                                               2, -2, 1, -2, 0, -2, 3, -2, /*  -7: (24)     p14 */
606b5a892a1SMatthew G. Knepley                                               1, -2, 0, -2, 2, -2, 3, -1, /*  -6: (34)     p6  */
607b5a892a1SMatthew G. Knepley                                               1, -1, 3, -3, 0, -3, 2, -2, /*  -5: (1243)   p10 */
608b5a892a1SMatthew G. Knepley                                               1, -3, 2, -1, 3, -1, 0, -3, /*  -4: (1432)   p9  */
609b5a892a1SMatthew G. Knepley                                               0, -3, 1, -1, 3, -3, 2, -3, /*  -3: (12)     p1  */
610b5a892a1SMatthew G. Knepley                                               0, -2, 2, -2, 1, -2, 3, -3, /*  -2: (23)     p2  */
611b5a892a1SMatthew G. Knepley                                               0, -1, 3, -2, 2, -3, 1, -2, /*  -1: (13)     p5  */
612b5a892a1SMatthew G. Knepley                                               0, 0,  1, 0,  2, 0,  3, 0,  /*   0: ()       p0  */
613b5a892a1SMatthew G. Knepley                                               0, 1,  3, 1,  1, 2,  2, 0,  /*   1: (123)    p4  */
614b5a892a1SMatthew G. Knepley                                               0, 2,  2, 1,  3, 0,  1, 2,  /*   2: (132)    p3  */
615b5a892a1SMatthew G. Knepley                                               1, 2,  0, 1,  3, 1,  2, 2,  /*   3: (12)(34) p7  */
616b5a892a1SMatthew G. Knepley                                               1, 0,  2, 0,  0, 0,  3, 1,  /*   4: (243)    p8  */
617b5a892a1SMatthew G. Knepley                                               1, 1,  3, 2,  2, 2,  0, 0,  /*   5: (143)    p11 */
618b5a892a1SMatthew G. Knepley                                               2, 1,  3, 0,  0, 2,  1, 0,  /*   6: (13)(24) p16 */
619b5a892a1SMatthew G. Knepley                                               2, 2,  1, 1,  3, 2,  0, 2,  /*   7: (142)    p15 */
620b5a892a1SMatthew G. Knepley                                               2, 0,  0, 0,  1, 0,  3, 2,  /*   8: (234)    p12 */
621b5a892a1SMatthew G. Knepley                                               3, 2,  2, 2,  1, 1,  0, 1,  /*   9: (14)(23) p23 */
622b5a892a1SMatthew G. Knepley                                               3, 0,  0, 2,  2, 1,  1, 1,  /*  10: (134)    p19 */
623b5a892a1SMatthew G. Knepley                                               3, 1,  1, 2,  0, 1,  2, 1 /*  11: (124)    p20 */};
624b5a892a1SMatthew G. Knepley   /* Each rotation determines a permutation of the four diagonals, and this defines the isomorphism with S_4 */
625ef8b56bfSJed Brown   static const PetscInt hexArr[48 * 6 * 2] = {
626b5a892a1SMatthew G. Knepley     2, -3, 3, -2, 4, -2, 5, -3, 1, -3, 0, -1, /* -24: reflect bottom and use -3 on top */
627b5a892a1SMatthew G. Knepley     4, -2, 5, -2, 0, -1, 1, -4, 3, -2, 2, -3, /* -23: reflect bottom and use -3 on top */
628b5a892a1SMatthew G. Knepley     5, -3, 4, -1, 1, -2, 0, -3, 3, -4, 2, -1, /* -22: reflect bottom and use -3 on top */
629b5a892a1SMatthew G. Knepley     3, -1, 2, -4, 4, -4, 5, -1, 0, -4, 1, -4, /* -21: reflect bottom and use -3 on top */
630b5a892a1SMatthew G. Knepley     3, -3, 2, -2, 5, -1, 4, -4, 1, -1, 0, -3, /* -20: reflect bottom and use -3 on top */
631b5a892a1SMatthew G. Knepley     4, -4, 5, -4, 1, -4, 0, -1, 2, -4, 3, -1, /* -19: reflect bottom and use -3 on top */
632b5a892a1SMatthew G. Knepley     2, -1, 3, -4, 5, -3, 4, -2, 0, -2, 1, -2, /* -18: reflect bottom and use -3 on top */
633b5a892a1SMatthew G. Knepley     5, -1, 4, -3, 0, -3, 1, -2, 2, -2, 3, -3, /* -17: reflect bottom and use -3 on top */
634b5a892a1SMatthew G. Knepley     4, -3, 5, -1, 3, -2, 2, -4, 1, -4, 0, -4, /* -16: reflect bottom and use -3 on top */
635b5a892a1SMatthew G. Knepley     5, -4, 4, -4, 3, -4, 2, -2, 0, -3, 1, -1, /* -15: reflect bottom and use -3 on top */
636b5a892a1SMatthew G. Knepley     3, -4, 2, -1, 1, -1, 0, -4, 4, -4, 5, -4, /* -14: reflect bottom and use -3 on top */
637b5a892a1SMatthew G. Knepley     2, -2, 3, -3, 0, -2, 1, -3, 4, -2, 5, -2, /* -13: reflect bottom and use -3 on top */
638b5a892a1SMatthew G. Knepley     1, -3, 0, -1, 4, -1, 5, -4, 3, -1, 2, -4, /* -12: reflect bottom and use -3 on top */
639b5a892a1SMatthew G. Knepley     1, -1, 0, -3, 5, -4, 4, -1, 2, -1, 3, -4, /* -11: reflect bottom and use -3 on top */
640b5a892a1SMatthew G. Knepley     5, -2, 4, -2, 2, -2, 3, -4, 1, -2, 0, -2, /* -10: reflect bottom and use -3 on top */
641b5a892a1SMatthew G. Knepley     1, -2, 0, -2, 2, -1, 3, -1, 4, -1, 5, -3, /*  -9: reflect bottom and use -3 on top */
642b5a892a1SMatthew G. Knepley     4, -1, 5, -3, 2, -4, 3, -2, 0, -1, 1, -3, /*  -8: reflect bottom and use -3 on top */
643b5a892a1SMatthew G. Knepley     3, -2, 2, -3, 0, -4, 1, -1, 5, -1, 4, -3, /*  -7: reflect bottom and use -3 on top */
644b5a892a1SMatthew G. Knepley     1, -4, 0, -4, 3, -1, 2, -1, 5, -4, 4, -4, /*  -6: reflect bottom and use -3 on top */
645b5a892a1SMatthew G. Knepley     2, -4, 3, -1, 1, -3, 0, -2, 5, -3, 4, -1, /*  -5: reflect bottom and use -3 on top */
646b5a892a1SMatthew G. Knepley     0, -4, 1, -4, 4, -3, 5, -2, 2, -3, 3, -2, /*  -4: reflect bottom and use -3 on top */
647b5a892a1SMatthew G. Knepley     0, -3, 1, -1, 3, -3, 2, -3, 4, -3, 5, -1, /*  -3: reflect bottom and use -3 on top */
648b5a892a1SMatthew G. Knepley     0, -2, 1, -2, 5, -2, 4, -3, 3, -3, 2, -2, /*  -2: reflect bottom and use -3 on top */
649b5a892a1SMatthew G. Knepley     0, -1, 1, -3, 2, -3, 3, -3, 5, -2, 4, -2, /*  -1: reflect bottom and use -3 on top */
650b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  5, 0,  /*   0: identity */
651b5a892a1SMatthew G. Knepley     0, 1,  1, 3,  5, 3,  4, 0,  2, 0,  3, 1,  /*   1: 90  rotation about z */
652b5a892a1SMatthew G. Knepley     0, 2,  1, 2,  3, 0,  2, 0,  5, 3,  4, 1,  /*   2: 180 rotation about z */
653b5a892a1SMatthew G. Knepley     0, 3,  1, 1,  4, 0,  5, 3,  3, 0,  2, 1,  /*   3: 270 rotation about z */
654b5a892a1SMatthew G. Knepley     2, 3,  3, 2,  1, 0,  0, 3,  4, 3,  5, 1,  /*   4: 90  rotation about x */
655b5a892a1SMatthew G. Knepley     1, 3,  0, 1,  3, 2,  2, 2,  4, 2,  5, 2,  /*   5: 180 rotation about x */
656b5a892a1SMatthew G. Knepley     3, 1,  2, 0,  0, 1,  1, 2,  4, 1,  5, 3,  /*   6: 270 rotation about x */
657b5a892a1SMatthew G. Knepley     4, 0,  5, 0,  2, 1,  3, 3,  1, 1,  0, 3,  /*   7: 90  rotation about y */
658b5a892a1SMatthew G. Knepley     1, 1,  0, 3,  2, 2,  3, 2,  5, 1,  4, 3,  /*   8: 180 rotation about y */
659b5a892a1SMatthew G. Knepley     5, 1,  4, 3,  2, 3,  3, 1,  0, 0,  1, 0,  /*   9: 270 rotation about y */
660b5a892a1SMatthew G. Knepley     1, 0,  0, 0,  5, 1,  4, 2,  3, 2,  2, 3,  /*  10: 180 rotation about x+y */
661b5a892a1SMatthew G. Knepley     1, 2,  0, 2,  4, 2,  5, 1,  2, 2,  3, 3,  /*  11: 180 rotation about x-y */
662b5a892a1SMatthew G. Knepley     2, 1,  3, 0,  0, 3,  1, 0,  5, 0,  4, 0,  /*  12: 180 rotation about y+z */
663b5a892a1SMatthew G. Knepley     3, 3,  2, 2,  1, 2,  0, 1,  5, 2,  4, 2,  /*  13: 180 rotation about y-z */
664b5a892a1SMatthew G. Knepley     5, 3,  4, 1,  3, 1,  2, 3,  1, 3,  0, 1,  /*  14: 180 rotation about z+x */
665b5a892a1SMatthew G. Knepley     4, 2,  5, 2,  3, 3,  2, 1,  0, 2,  1, 2,  /*  15: 180 rotation about z-x */
666b5a892a1SMatthew G. Knepley     5, 0,  4, 0,  0, 0,  1, 3,  3, 1,  2, 0,  /*  16: 120 rotation about x+y+z (v0v6) */
667b5a892a1SMatthew G. Knepley     2, 0,  3, 1,  5, 0,  4, 3,  1, 0,  0, 0,  /*  17: 240 rotation about x+y+z (v0v6) */
668b5a892a1SMatthew G. Knepley     4, 3,  5, 1,  1, 1,  0, 2,  3, 3,  2, 2,  /*  18: 120 rotation about x+y-z (v4v2) */
669b5a892a1SMatthew G. Knepley     3, 2,  2, 3,  5, 2,  4, 1,  0, 1,  1, 3,  /*  19: 240 rotation about x+y-z (v4v2) */
670b5a892a1SMatthew G. Knepley     3, 0,  2, 1,  4, 1,  5, 2,  1, 2,  0, 2,  /*  20: 120 rotation about x-y+z (v1v5) */
671b5a892a1SMatthew G. Knepley     5, 2,  4, 2,  1, 3,  0, 0,  2, 3,  3, 2,  /*  21: 240 rotation about x-y+z (v1v5) */
672b5a892a1SMatthew G. Knepley     4, 1,  5, 3,  0, 2,  1, 1,  2, 1,  3, 0,  /*  22: 120 rotation about x-y-z (v7v3) */
673b5a892a1SMatthew G. Knepley     2, 2,  3, 3,  4, 3,  5, 0,  0, 3,  1, 1,  /*  23: 240 rotation about x-y-z (v7v3) */
674b5a892a1SMatthew G. Knepley   };
675ef8b56bfSJed Brown   static const PetscInt tripArr[12 * 5 * 2] = {
676b5a892a1SMatthew G. Knepley     1, -3, 0, -1, 3, -1, 4, -1, 2, -1, /* -6: reflect bottom and top */
677b5a892a1SMatthew G. Knepley     1, -1, 0, -3, 4, -1, 2, -1, 3, -1, /* -5: reflect bottom and top */
678b5a892a1SMatthew G. Knepley     1, -2, 0, -2, 2, -1, 3, -1, 4, -1, /* -4: reflect bottom and top */
679b5a892a1SMatthew G. Knepley     0, -3, 1, -1, 3, -3, 2, -3, 4, -3, /* -3: reflect bottom and top */
680b5a892a1SMatthew G. Knepley     0, -2, 1, -2, 4, -3, 3, -3, 2, -3, /* -2: reflect bottom and top */
681b5a892a1SMatthew G. Knepley     0, -1, 1, -3, 2, -3, 4, -3, 3, -3, /* -1: reflect bottom and top */
682b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  /*  0: identity */
683b5a892a1SMatthew G. Knepley     0, 1,  1, 2,  4, 0,  2, 0,  3, 0,  /*  1: 120 rotation about z */
684b5a892a1SMatthew G. Knepley     0, 2,  1, 1,  3, 0,  4, 0,  2, 0,  /*  2: 240 rotation about z */
685b5a892a1SMatthew G. Knepley     1, 1,  0, 2,  2, 2,  4, 2,  3, 2,  /*  3: 180 rotation about y of 0 */
686b5a892a1SMatthew G. Knepley     1, 0,  0, 0,  4, 2,  3, 2,  2, 2,  /*  4: 180 rotation about y of 1 */
687b5a892a1SMatthew G. Knepley     1, 2,  0, 1,  3, 2,  2, 2,  4, 2,  /*  5: 180 rotation about y of 2 */
688b5a892a1SMatthew G. Knepley   };
689b5a892a1SMatthew G. Knepley   /* a: rotate 120 about z
690b5a892a1SMatthew G. Knepley      b: swap top and bottom segments
691b5a892a1SMatthew G. Knepley      r: reflect */
692ef8b56bfSJed Brown   static const PetscInt ttriArr[12 * 5 * 2] = {
693b5a892a1SMatthew G. Knepley     1, -3, 0, -3, 2, -2, 4, -2, 3, -2, /* -6: r b a^2 */
694b5a892a1SMatthew G. Knepley     1, -2, 0, -2, 4, -2, 3, -2, 2, -2, /* -5: r b a */
695b5a892a1SMatthew G. Knepley     1, -1, 0, -1, 3, -2, 2, -2, 4, -2, /* -4: r b */
696b5a892a1SMatthew G. Knepley     0, -3, 1, -3, 2, -1, 4, -1, 3, -1, /* -3: r a^2 */
697b5a892a1SMatthew G. Knepley     0, -2, 1, -2, 4, -1, 3, -1, 2, -1, /* -2: r a */
698b5a892a1SMatthew G. Knepley     0, -1, 1, -1, 3, -1, 2, -1, 4, -1, /* -1: r */
699b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  /*  0: identity */
700b5a892a1SMatthew G. Knepley     0, 1,  1, 1,  3, 0,  4, 0,  2, 0,  /*  1: a */
701b5a892a1SMatthew G. Knepley     0, 2,  1, 2,  4, 0,  2, 0,  3, 0,  /*  2: a^2 */
702b5a892a1SMatthew G. Knepley     1, 0,  0, 0,  2, 1,  3, 1,  4, 1,  /*  3: b */
703b5a892a1SMatthew G. Knepley     1, 1,  0, 1,  3, 1,  4, 1,  2, 1,  /*  4: b a */
704b5a892a1SMatthew G. Knepley     1, 2,  0, 2,  4, 1,  2, 1,  3, 1,  /*  5: b a^2 */
705b5a892a1SMatthew G. Knepley   };
706b5a892a1SMatthew G. Knepley   /* a: rotate 90 about z
707b5a892a1SMatthew G. Knepley      b: swap top and bottom segments
708b5a892a1SMatthew G. Knepley      r: reflect */
709ef8b56bfSJed Brown   static const PetscInt tquadArr[16 * 6 * 2] = {
710b5a892a1SMatthew G. Knepley     1, -4, 0, -4, 3, -2, 2, -2, 5, -2, 4, -2, /* -8: r b a^3 */
711b5a892a1SMatthew G. Knepley     1, -3, 0, -3, 2, -2, 5, -2, 4, -2, 3, -2, /* -7: r b a^2 */
712b5a892a1SMatthew G. Knepley     1, -2, 0, -2, 5, -2, 4, -2, 3, -2, 2, -2, /* -6: r b a */
713b5a892a1SMatthew G. Knepley     1, -1, 0, -1, 4, -2, 3, -2, 2, -2, 5, -2, /* -5: r b */
714b5a892a1SMatthew G. Knepley     0, -4, 1, -4, 3, -1, 2, -1, 5, -1, 4, -1, /* -4: r a^3 */
715b5a892a1SMatthew G. Knepley     0, -3, 1, -3, 2, -1, 5, -1, 4, -1, 3, -1, /* -3: r a^2 */
716b5a892a1SMatthew G. Knepley     0, -2, 1, -2, 5, -1, 4, -1, 3, -1, 2, -1, /* -2: r a */
717b5a892a1SMatthew G. Knepley     0, -1, 1, -1, 4, -1, 3, -1, 2, -1, 5, -1, /* -1: r */
718b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  5, 0,  /*  0: identity */
719b5a892a1SMatthew G. Knepley     0, 1,  1, 1,  3, 0,  4, 0,  5, 0,  2, 0,  /*  1: a */
720b5a892a1SMatthew G. Knepley     0, 2,  1, 2,  4, 0,  5, 0,  2, 0,  3, 0,  /*  2: a^2 */
721b5a892a1SMatthew G. Knepley     0, 3,  1, 3,  5, 0,  2, 0,  3, 0,  4, 0,  /*  3: a^3 */
722b5a892a1SMatthew G. Knepley     1, 0,  0, 0,  2, 1,  3, 1,  4, 1,  5, 1,  /*  4: b */
723b5a892a1SMatthew G. Knepley     1, 1,  0, 1,  3, 1,  4, 1,  5, 1,  2, 1,  /*  5: b a */
724b5a892a1SMatthew G. Knepley     1, 2,  0, 2,  4, 1,  5, 1,  2, 1,  3, 1,  /*  6: b a^2 */
725b5a892a1SMatthew G. Knepley     1, 3,  0, 3,  5, 1,  2, 1,  3, 1,  4, 1,  /*  7: b a^3 */
726b5a892a1SMatthew G. Knepley   };
727ef8b56bfSJed Brown   static const PetscInt pyrArr[8 * 5 * 2] = {
728b5a892a1SMatthew G. Knepley     0, -4, 2, -3, 1, -3, 4, -3, 3, -3, /* -4: Reflect bottom face */
729b5a892a1SMatthew G. Knepley     0, -3, 3, -3, 2, -3, 1, -3, 4, -3, /* -3: Reflect bottom face */
730b5a892a1SMatthew G. Knepley     0, -2, 4, -3, 3, -3, 2, -3, 1, -3, /* -2: Reflect bottom face */
731b5a892a1SMatthew G. Knepley     0, -1, 1, -3, 4, -3, 3, -3, 2, -3, /* -1: Reflect bottom face */
732b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  /*  0: identity */
733b5a892a1SMatthew G. Knepley     0, 1,  4, 0,  1, 0,  2, 0,  3, 0,  /*  1:  90 rotation about z */
734b5a892a1SMatthew G. Knepley     0, 2,  3, 0,  4, 0,  1, 0,  2, 0,  /*  2: 180 rotation about z */
735b5a892a1SMatthew G. Knepley     0, 3,  2, 0,  3, 0,  4, 0,  1, 0,  /*  3: 270 rotation about z */
736b5a892a1SMatthew G. Knepley   };
737b5a892a1SMatthew G. Knepley   switch (ct) {
738d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
739d71ae5a4SJacob Faibussowitsch     return pntArr;
740d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
741d71ae5a4SJacob Faibussowitsch     return &segArr[(o + 1) * 2 * 2];
742d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
743d71ae5a4SJacob Faibussowitsch     return &segArr[(o + 1) * 2 * 2];
744d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
745d71ae5a4SJacob Faibussowitsch     return &triArr[(o + 3) * 3 * 2];
746d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
747d71ae5a4SJacob Faibussowitsch     return &quadArr[(o + 4) * 4 * 2];
748d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
749d71ae5a4SJacob Faibussowitsch     return &tsegArr[(o + 2) * 4 * 2];
750d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
751d71ae5a4SJacob Faibussowitsch     return &tetArr[(o + 12) * 4 * 2];
752d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
753d71ae5a4SJacob Faibussowitsch     return &hexArr[(o + 24) * 6 * 2];
754d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
755d71ae5a4SJacob Faibussowitsch     return &tripArr[(o + 6) * 5 * 2];
756d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
757d71ae5a4SJacob Faibussowitsch     return &ttriArr[(o + 6) * 5 * 2];
758d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
759d71ae5a4SJacob Faibussowitsch     return &tquadArr[(o + 8) * 6 * 2];
760d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
761d71ae5a4SJacob Faibussowitsch     return &pyrArr[(o + 4) * 5 * 2];
762d71ae5a4SJacob Faibussowitsch   default:
763f22e26b7SPierre Jolivet     return PETSC_NULLPTR;
764b5a892a1SMatthew G. Knepley   }
765b5a892a1SMatthew G. Knepley }
766b5a892a1SMatthew G. Knepley 
767d5b43468SJose E. Roman /* A vertex arrangement is a vertex order */
768d71ae5a4SJacob Faibussowitsch static inline const PetscInt *DMPolytopeTypeGetVertexArrangment(DMPolytopeType ct, PetscInt o)
769d71ae5a4SJacob Faibussowitsch {
770ef8b56bfSJed Brown   static const PetscInt pntVerts[1]      = {0};
7719371c9d4SSatish Balay   static const PetscInt segVerts[2 * 2]  = {1, 0, 0, 1};
7729371c9d4SSatish Balay   static const PetscInt triVerts[6 * 3]  = {1, 0, 2, 0, 2, 1, 2, 1, 0, 0, 1, 2, 1, 2, 0, 2, 0, 1};
7739371c9d4SSatish Balay   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};
7749371c9d4SSatish Balay   static const PetscInt tsegVerts[4 * 4] = {3, 2, 1, 0, 1, 0, 3, 2, 0, 1, 2, 3, 2, 3, 0, 1};
7759371c9d4SSatish Balay   static const PetscInt tetVerts[24 * 4] = {2, 3, 1, 0, /* -12: (1324)   p22 */
776b5a892a1SMatthew G. Knepley                                             3, 1, 2, 0, /* -11: (14)     p21 */
777b5a892a1SMatthew G. Knepley                                             1, 2, 3, 0, /* -10: (1234)   p18 */
778b5a892a1SMatthew G. Knepley                                             3, 2, 0, 1, /*  -9: (1423)   p17 */
779b5a892a1SMatthew G. Knepley                                             2, 0, 3, 1, /*  -8: (1342)   p13 */
780b5a892a1SMatthew G. Knepley                                             0, 3, 2, 1, /*  -7: (24)     p14 */
781b5a892a1SMatthew G. Knepley                                             0, 1, 3, 2, /*  -6: (34)     p6  */
782b5a892a1SMatthew G. Knepley                                             1, 3, 0, 2, /*  -5: (1243)   p10 */
783b5a892a1SMatthew G. Knepley                                             3, 0, 1, 2, /*  -4: (1432    p9  */
784b5a892a1SMatthew G. Knepley                                             1, 0, 2, 3, /*  -3: (12)     p1  */
785b5a892a1SMatthew G. Knepley                                             0, 2, 1, 3, /*  -2: (23)     p2  */
786b5a892a1SMatthew G. Knepley                                             2, 1, 0, 3, /*  -1: (13)     p5  */
787b5a892a1SMatthew G. Knepley                                             0, 1, 2, 3, /*   0: ()       p0  */
788b5a892a1SMatthew G. Knepley                                             1, 2, 0, 3, /*   1: (123)    p4  */
789b5a892a1SMatthew G. Knepley                                             2, 0, 1, 3, /*   2: (132)    p3  */
790b5a892a1SMatthew G. Knepley                                             1, 0, 3, 2, /*   3: (12)(34) p7  */
791b5a892a1SMatthew G. Knepley                                             0, 3, 1, 2, /*   4: (243)    p8  */
792b5a892a1SMatthew G. Knepley                                             3, 1, 0, 2, /*   5: (143)    p11 */
793b5a892a1SMatthew G. Knepley                                             2, 3, 0, 1, /*   6: (13)(24) p16 */
794b5a892a1SMatthew G. Knepley                                             3, 0, 2, 1, /*   7: (142)    p15 */
795b5a892a1SMatthew G. Knepley                                             0, 2, 3, 1, /*   8: (234)    p12 */
796b5a892a1SMatthew G. Knepley                                             3, 2, 1, 0, /*   9: (14)(23) p23 */
797b5a892a1SMatthew G. Knepley                                             2, 1, 3, 0, /*  10: (134)    p19 */
798b5a892a1SMatthew G. Knepley                                             1, 3, 2, 0 /*  11: (124)    p20 */};
799ef8b56bfSJed Brown   static const PetscInt hexVerts[48 * 8] = {
800b5a892a1SMatthew G. Knepley     3, 0, 4, 5, 2, 6, 7, 1, /* -24: reflected 23 */
801b5a892a1SMatthew G. Knepley     3, 5, 6, 2, 0, 1, 7, 4, /* -23: reflected 22 */
802b5a892a1SMatthew G. Knepley     4, 0, 1, 7, 5, 6, 2, 3, /* -22: reflected 21 */
803b5a892a1SMatthew G. Knepley     6, 7, 1, 2, 5, 3, 0, 4, /* -21: reflected 20 */
804b5a892a1SMatthew G. Knepley     1, 2, 6, 7, 0, 4, 5, 3, /* -20: reflected 19 */
805b5a892a1SMatthew G. Knepley     6, 2, 3, 5, 7, 4, 0, 1, /* -19: reflected 18 */
806b5a892a1SMatthew G. Knepley     4, 5, 3, 0, 7, 1, 2, 6, /* -18: reflected 17 */
807b5a892a1SMatthew G. Knepley     1, 7, 4, 0, 2, 3, 5, 6, /* -17: reflected 16 */
808b5a892a1SMatthew G. Knepley     2, 3, 5, 6, 1, 7, 4, 0, /* -16: reflected 15 */
809b5a892a1SMatthew G. Knepley     7, 4, 0, 1, 6, 2, 3, 5, /* -15: reflected 14 */
810b5a892a1SMatthew G. Knepley     7, 1, 2, 6, 4, 5, 3, 0, /* -14: reflected 13 */
811b5a892a1SMatthew G. Knepley     0, 4, 5, 3, 1, 2, 6, 7, /* -13: reflected 12 */
812b5a892a1SMatthew G. Knepley     5, 4, 7, 6, 3, 2, 1, 0, /* -12: reflected 11 */
813b5a892a1SMatthew G. Knepley     7, 6, 5, 4, 1, 0, 3, 2, /* -11: reflected 10 */
814b5a892a1SMatthew G. Knepley     0, 1, 7, 4, 3, 5, 6, 2, /* -10: reflected  9 */
815b5a892a1SMatthew G. Knepley     4, 7, 6, 5, 0, 3, 2, 1, /*  -9: reflected  8 */
816b5a892a1SMatthew G. Knepley     5, 6, 2, 3, 4, 0, 1, 7, /*  -8: reflected  7 */
817b5a892a1SMatthew G. Knepley     2, 6, 7, 1, 3, 0, 4, 5, /*  -7: reflected  6 */
818b5a892a1SMatthew G. Knepley     6, 5, 4, 7, 2, 1, 0, 3, /*  -6: reflected  5 */
819b5a892a1SMatthew G. Knepley     5, 3, 0, 4, 6, 7, 1, 2, /*  -5: reflected  4 */
820b5a892a1SMatthew G. Knepley     2, 1, 0, 3, 6, 5, 4, 7, /*  -4: reflected  3 */
821b5a892a1SMatthew G. Knepley     1, 0, 3, 2, 7, 6, 5, 4, /*  -3: reflected  2 */
822b5a892a1SMatthew G. Knepley     0, 3, 2, 1, 4, 7, 6, 5, /*  -2: reflected  1 */
823b5a892a1SMatthew G. Knepley     3, 2, 1, 0, 5, 4, 7, 6, /*  -1: reflected  0 */
824b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, 5, 6, 7, /*   0: identity */
825b5a892a1SMatthew G. Knepley     1, 2, 3, 0, 7, 4, 5, 6, /*   1: 90  rotation about z */
826b5a892a1SMatthew G. Knepley     2, 3, 0, 1, 6, 7, 4, 5, /*   2: 180 rotation about z */
827b5a892a1SMatthew G. Knepley     3, 0, 1, 2, 5, 6, 7, 4, /*   3: 270 rotation about z */
828b5a892a1SMatthew G. Knepley     4, 0, 3, 5, 7, 6, 2, 1, /*   4: 90  rotation about x */
829b5a892a1SMatthew G. Knepley     7, 4, 5, 6, 1, 2, 3, 0, /*   5: 180 rotation about x */
830b5a892a1SMatthew G. Knepley     1, 7, 6, 2, 0, 3, 5, 4, /*   6: 270 rotation about x */
831b5a892a1SMatthew G. Knepley     3, 2, 6, 5, 0, 4, 7, 1, /*   7: 90  rotation about y */
832b5a892a1SMatthew G. Knepley     5, 6, 7, 4, 3, 0, 1, 2, /*   8: 180 rotation about y */
833b5a892a1SMatthew G. Knepley     4, 7, 1, 0, 5, 3, 2, 6, /*   9: 270 rotation about y */
834b5a892a1SMatthew G. Knepley     4, 5, 6, 7, 0, 1, 2, 3, /*  10: 180 rotation about x+y */
835b5a892a1SMatthew G. Knepley     6, 7, 4, 5, 2, 3, 0, 1, /*  11: 180 rotation about x-y */
836b5a892a1SMatthew G. Knepley     3, 5, 4, 0, 2, 1, 7, 6, /*  12: 180 rotation about y+z */
837b5a892a1SMatthew G. Knepley     6, 2, 1, 7, 5, 4, 0, 3, /*  13: 180 rotation about y-z */
838b5a892a1SMatthew G. Knepley     1, 0, 4, 7, 2, 6, 5, 3, /*  14: 180 rotation about z+x */
839b5a892a1SMatthew G. Knepley     6, 5, 3, 2, 7, 1, 0, 4, /*  15: 180 rotation about z-x */
840b5a892a1SMatthew G. Knepley     0, 4, 7, 1, 3, 2, 6, 5, /*  16: 120 rotation about x+y+z (v0v6) */
841b5a892a1SMatthew G. Knepley     0, 3, 5, 4, 1, 7, 6, 2, /*  17: 240 rotation about x+y+z (v0v6) */
842b5a892a1SMatthew G. Knepley     5, 3, 2, 6, 4, 7, 1, 0, /*  18: 120 rotation about x+y-z (v4v2) */
843b5a892a1SMatthew G. Knepley     7, 6, 2, 1, 4, 0, 3, 5, /*  19: 240 rotation about x+y-z (v4v2) */
844b5a892a1SMatthew G. Knepley     2, 1, 7, 6, 3, 5, 4, 0, /*  20: 120 rotation about x-y+z (v1v5) */
845b5a892a1SMatthew G. Knepley     7, 1, 0, 4, 6, 5, 3, 2, /*  21: 240 rotation about x-y+z (v1v5) */
846b5a892a1SMatthew G. Knepley     2, 6, 5, 3, 1, 0, 4, 7, /*  22: 120 rotation about x-y-z (v7v3) */
847b5a892a1SMatthew G. Knepley     5, 4, 0, 3, 6, 2, 1, 7, /*  23: 240 rotation about x-y-z (v7v3) */
848b5a892a1SMatthew G. Knepley   };
849ef8b56bfSJed Brown   static const PetscInt tripVerts[12 * 6] = {
850b5a892a1SMatthew G. Knepley     4, 3, 5, 2, 1, 0, /* -6: reflect bottom and top */
851b5a892a1SMatthew G. Knepley     5, 4, 3, 1, 0, 2, /* -5: reflect bottom and top */
852b5a892a1SMatthew G. Knepley     3, 5, 4, 0, 2, 1, /* -4: reflect bottom and top */
853b5a892a1SMatthew G. Knepley     1, 0, 2, 5, 4, 3, /* -3: reflect bottom and top */
854b5a892a1SMatthew G. Knepley     0, 2, 1, 3, 5, 4, /* -2: reflect bottom and top */
855b5a892a1SMatthew G. Knepley     2, 1, 0, 4, 3, 5, /* -1: reflect bottom and top */
856b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, 5, /*  0: identity */
857b5a892a1SMatthew G. Knepley     1, 2, 0, 5, 3, 4, /*  1: 120 rotation about z */
858b5a892a1SMatthew G. Knepley     2, 0, 1, 4, 5, 3, /*  2: 240 rotation about z */
859b5a892a1SMatthew G. Knepley     4, 5, 3, 2, 0, 1, /*  3: 180 rotation about y of 0 */
860b5a892a1SMatthew G. Knepley     3, 4, 5, 0, 1, 2, /*  4: 180 rotation about y of 1 */
861b5a892a1SMatthew G. Knepley     5, 3, 4, 1, 2, 0, /*  5: 180 rotation about y of 2 */
862b5a892a1SMatthew G. Knepley   };
863ef8b56bfSJed Brown   static const PetscInt ttriVerts[12 * 6] = {
864b5a892a1SMatthew G. Knepley     4, 3, 5, 1, 0, 2, /* -6: r b a^2 */
865b5a892a1SMatthew G. Knepley     3, 5, 4, 0, 2, 1, /* -5: r b a */
866b5a892a1SMatthew G. Knepley     5, 4, 3, 2, 1, 0, /* -4: r b */
867b5a892a1SMatthew G. Knepley     1, 0, 2, 4, 3, 5, /* -3: r a^2 */
868b5a892a1SMatthew G. Knepley     0, 2, 1, 3, 5, 4, /* -2: r a */
869b5a892a1SMatthew G. Knepley     2, 1, 0, 5, 4, 3, /* -1: r */
870b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, 5, /*  0: identity */
871b5a892a1SMatthew G. Knepley     1, 2, 0, 4, 5, 3, /*  1: a */
872b5a892a1SMatthew G. Knepley     2, 0, 1, 5, 3, 4, /*  2: a^2 */
873b5a892a1SMatthew G. Knepley     3, 4, 5, 0, 1, 2, /*  3: b */
874b5a892a1SMatthew G. Knepley     4, 5, 3, 1, 2, 0, /*  4: b a */
875b5a892a1SMatthew G. Knepley     5, 3, 4, 2, 0, 1, /*  5: b a^2 */
876b5a892a1SMatthew G. Knepley   };
877b5a892a1SMatthew G. Knepley   /* a: rotate 90 about z
878b5a892a1SMatthew G. Knepley      b: swap top and bottom segments
879b5a892a1SMatthew G. Knepley      r: reflect */
880ef8b56bfSJed Brown   static const PetscInt tquadVerts[16 * 8] = {
881b5a892a1SMatthew G. Knepley     6, 5, 4, 7, 2, 1, 0, 3, /* -8: r b a^3 */
882b5a892a1SMatthew G. Knepley     5, 4, 7, 6, 1, 0, 3, 2, /* -7: r b a^2 */
883b5a892a1SMatthew G. Knepley     4, 7, 6, 5, 0, 3, 2, 1, /* -6: r b a */
884b5a892a1SMatthew G. Knepley     7, 6, 5, 4, 3, 2, 1, 0, /* -5: r b */
885b5a892a1SMatthew G. Knepley     2, 1, 0, 3, 6, 5, 4, 7, /* -4: r a^3 */
886b5a892a1SMatthew G. Knepley     1, 0, 3, 2, 5, 4, 7, 6, /* -3: r a^2 */
887b5a892a1SMatthew G. Knepley     0, 3, 2, 1, 4, 7, 6, 5, /* -2: r a */
888b5a892a1SMatthew G. Knepley     3, 2, 1, 0, 7, 6, 5, 4, /* -1: r */
889b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, 5, 6, 7, /*  0: identity */
890b5a892a1SMatthew G. Knepley     1, 2, 3, 0, 5, 6, 7, 4, /*  1: a */
891b5a892a1SMatthew G. Knepley     2, 3, 0, 1, 6, 7, 4, 5, /*  2: a^2 */
892b5a892a1SMatthew G. Knepley     3, 0, 1, 2, 7, 4, 5, 6, /*  3: a^3 */
893b5a892a1SMatthew G. Knepley     4, 5, 6, 7, 0, 1, 2, 3, /*  4: b */
894b5a892a1SMatthew G. Knepley     5, 6, 7, 4, 1, 2, 3, 0, /*  5: b a */
895b5a892a1SMatthew G. Knepley     6, 7, 4, 5, 2, 3, 0, 1, /*  6: b a^2 */
896b5a892a1SMatthew G. Knepley     7, 4, 5, 6, 3, 0, 1, 2, /*  7: b a^3 */
897b5a892a1SMatthew G. Knepley   };
898ef8b56bfSJed Brown   static const PetscInt pyrVerts[8 * 5] = {
899b5a892a1SMatthew G. Knepley     2, 1, 0, 3, 4, /* -4: Reflect bottom face */
900b5a892a1SMatthew G. Knepley     1, 0, 3, 2, 4, /* -3: Reflect bottom face */
901b5a892a1SMatthew G. Knepley     0, 3, 2, 1, 4, /* -2: Reflect bottom face */
902b5a892a1SMatthew G. Knepley     3, 2, 1, 0, 4, /* -1: Reflect bottom face */
903b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, /*  0: identity */
904b5a892a1SMatthew G. Knepley     1, 2, 3, 0, 4, /*  1:  90 rotation about z */
905b5a892a1SMatthew G. Knepley     2, 3, 0, 1, 4, /*  2: 180 rotation about z */
906b5a892a1SMatthew G. Knepley     3, 0, 1, 2, 4, /*  3: 270 rotation about z */
907b5a892a1SMatthew G. Knepley   };
908b5a892a1SMatthew G. Knepley   switch (ct) {
909d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
910d71ae5a4SJacob Faibussowitsch     return pntVerts;
911d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
912d71ae5a4SJacob Faibussowitsch     return &segVerts[(o + 1) * 2];
913d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
914d71ae5a4SJacob Faibussowitsch     return &segVerts[(o + 1) * 2];
915d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
916d71ae5a4SJacob Faibussowitsch     return &triVerts[(o + 3) * 3];
917d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
918d71ae5a4SJacob Faibussowitsch     return &quadVerts[(o + 4) * 4];
919d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
920d71ae5a4SJacob Faibussowitsch     return &tsegVerts[(o + 2) * 4];
921d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
922d71ae5a4SJacob Faibussowitsch     return &tetVerts[(o + 12) * 4];
923d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
924d71ae5a4SJacob Faibussowitsch     return &hexVerts[(o + 24) * 8];
925d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
926d71ae5a4SJacob Faibussowitsch     return &tripVerts[(o + 6) * 6];
927d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
928d71ae5a4SJacob Faibussowitsch     return &ttriVerts[(o + 6) * 6];
929d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
930d71ae5a4SJacob Faibussowitsch     return &tquadVerts[(o + 8) * 8];
931d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
932d71ae5a4SJacob Faibussowitsch     return &pyrVerts[(o + 4) * 5];
933d71ae5a4SJacob Faibussowitsch   default:
934f22e26b7SPierre Jolivet     return PETSC_NULLPTR;
935b5a892a1SMatthew G. Knepley   }
936b5a892a1SMatthew G. Knepley }
937b5a892a1SMatthew G. Knepley 
938b5a892a1SMatthew G. Knepley /* This is orientation o1 acting on orientation o2 */
939d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeComposeOrientation(DMPolytopeType ct, PetscInt o1, PetscInt o2)
940d71ae5a4SJacob Faibussowitsch {
9419371c9d4SSatish Balay   static const PetscInt segMult[2 * 2]   = {0, -1, -1, 0};
9429371c9d4SSatish Balay   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};
9439371c9d4SSatish Balay   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,
9449371c9d4SSatish Balay                                             -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};
9459371c9d4SSatish Balay   static const PetscInt tsegMult[4 * 4]  = {0, 1, -2, -1, 1, 0, -1, -2, -2, -1, 0, 1, -1, -2, 1, 0};
946ef8b56bfSJed Brown   static const PetscInt tetMult[24 * 24] = {
9479371c9d4SSatish Balay     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,
9489371c9d4SSatish Balay     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,
9499371c9d4SSatish Balay     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,
9509371c9d4SSatish Balay     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,
9519371c9d4SSatish Balay     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,
9529371c9d4SSatish Balay     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,
9539371c9d4SSatish Balay     -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,
9549371c9d4SSatish Balay     -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,
9559371c9d4SSatish Balay     -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,
9569371c9d4SSatish Balay     -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,
9579371c9d4SSatish Balay     -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,
9589371c9d4SSatish Balay     -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,
959b5a892a1SMatthew G. Knepley   };
960ef8b56bfSJed Brown   static const PetscInt hexMult[48 * 48] = {
961b5a892a1SMatthew G. Knepley     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,
962b5a892a1SMatthew G. Knepley     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,
963b5a892a1SMatthew G. Knepley     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,
964b5a892a1SMatthew G. Knepley     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,
965b5a892a1SMatthew G. Knepley     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,
966b5a892a1SMatthew G. Knepley     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,
967b5a892a1SMatthew G. Knepley     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,
968b5a892a1SMatthew G. Knepley     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,
969b5a892a1SMatthew G. Knepley     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,
970b5a892a1SMatthew G. Knepley     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,
971b5a892a1SMatthew G. Knepley     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,
972b5a892a1SMatthew G. Knepley     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,
973b5a892a1SMatthew G. Knepley     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,
974b5a892a1SMatthew G. Knepley     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,
975b5a892a1SMatthew G. Knepley     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,
976b5a892a1SMatthew G. Knepley     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,
977b5a892a1SMatthew G. Knepley     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,
978b5a892a1SMatthew G. Knepley     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,
979b5a892a1SMatthew G. Knepley     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,
980b5a892a1SMatthew G. Knepley     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,
981b5a892a1SMatthew G. Knepley     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,
982b5a892a1SMatthew G. Knepley     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,
983b5a892a1SMatthew G. Knepley     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,
984b5a892a1SMatthew G. Knepley     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,
985b5a892a1SMatthew G. Knepley     -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,
986b5a892a1SMatthew G. Knepley     -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,
987b5a892a1SMatthew G. Knepley     -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,
988b5a892a1SMatthew G. Knepley     -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,
989b5a892a1SMatthew G. Knepley     -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,
990b5a892a1SMatthew G. Knepley     -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,
991b5a892a1SMatthew G. Knepley     -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,
992b5a892a1SMatthew G. Knepley     -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,
993b5a892a1SMatthew G. Knepley     -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,
994b5a892a1SMatthew G. Knepley     -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,
995b5a892a1SMatthew G. Knepley     -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,
996b5a892a1SMatthew G. Knepley     -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,
997b5a892a1SMatthew G. Knepley     -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,
998b5a892a1SMatthew G. Knepley     -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,
999b5a892a1SMatthew G. Knepley     -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,
1000b5a892a1SMatthew G. Knepley     -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,
1001b5a892a1SMatthew G. Knepley     -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,
1002b5a892a1SMatthew G. Knepley     -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,
1003b5a892a1SMatthew G. Knepley     -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,
1004b5a892a1SMatthew G. Knepley     -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,
1005b5a892a1SMatthew G. Knepley     -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,
1006b5a892a1SMatthew G. Knepley     -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,
1007b5a892a1SMatthew G. Knepley     -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,
1008b5a892a1SMatthew G. Knepley     -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,
1009b5a892a1SMatthew G. Knepley   };
1010ef8b56bfSJed Brown   static const PetscInt tripMult[12 * 12] = {
10119371c9d4SSatish Balay     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,
10129371c9d4SSatish Balay     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,
10139371c9d4SSatish Balay     -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,
1014b5a892a1SMatthew G. Knepley   };
1015ef8b56bfSJed Brown   static const PetscInt ttriMult[12 * 12] = {
10169371c9d4SSatish Balay     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,
10179371c9d4SSatish Balay     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,
10189371c9d4SSatish Balay     -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,
1019b5a892a1SMatthew G. Knepley   };
1020ef8b56bfSJed Brown   static const PetscInt tquadMult[16 * 16] = {
10219371c9d4SSatish Balay     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,
10229371c9d4SSatish Balay     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,
10239371c9d4SSatish Balay     -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,
10249371c9d4SSatish Balay     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,
10259371c9d4SSatish Balay     -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,
1026b5a892a1SMatthew G. Knepley   };
1027ef8b56bfSJed Brown   static const PetscInt pyrMult[8 * 8] = {
10289371c9d4SSatish Balay     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,
1029b5a892a1SMatthew G. Knepley   };
1030b5a892a1SMatthew G. Knepley   switch (ct) {
1031d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
1032d71ae5a4SJacob Faibussowitsch     return 0;
1033b5a892a1SMatthew G. Knepley   case DM_POLYTOPE_SEGMENT:
1034d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
1035d71ae5a4SJacob Faibussowitsch     return segMult[(o1 + 1) * 2 + o2 + 1];
1036d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
1037d71ae5a4SJacob Faibussowitsch     return triMult[(o1 + 3) * 6 + o2 + 3];
1038d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
1039d71ae5a4SJacob Faibussowitsch     return quadMult[(o1 + 4) * 8 + o2 + 4];
1040d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
1041d71ae5a4SJacob Faibussowitsch     return tsegMult[(o1 + 2) * 4 + o2 + 2];
1042d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
1043d71ae5a4SJacob Faibussowitsch     return tetMult[(o1 + 12) * 24 + o2 + 12];
1044d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
1045d71ae5a4SJacob Faibussowitsch     return hexMult[(o1 + 24) * 48 + o2 + 24];
1046d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
1047d71ae5a4SJacob Faibussowitsch     return tripMult[(o1 + 6) * 12 + o2 + 6];
1048d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
1049d71ae5a4SJacob Faibussowitsch     return ttriMult[(o1 + 6) * 12 + o2 + 6];
1050d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
1051d71ae5a4SJacob Faibussowitsch     return tquadMult[(o1 + 8) * 16 + o2 + 8];
1052d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
1053d71ae5a4SJacob Faibussowitsch     return pyrMult[(o1 + 4) * 8 + o2 + 4];
1054d71ae5a4SJacob Faibussowitsch   default:
1055d71ae5a4SJacob Faibussowitsch     return 0;
1056b5a892a1SMatthew G. Knepley   }
1057b5a892a1SMatthew G. Knepley }
1058b5a892a1SMatthew G. Knepley 
1059b5a892a1SMatthew G. Knepley /* This is orientation o1 acting on orientation o2^{-1} */
1060d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeComposeOrientationInv(DMPolytopeType ct, PetscInt o1, PetscInt o2)
1061d71ae5a4SJacob Faibussowitsch {
1062ef8b56bfSJed Brown   static const PetscInt triInv[6]    = {-3, -2, -1, 0, 2, 1};
1063ef8b56bfSJed Brown   static const PetscInt quadInv[8]   = {-4, -3, -2, -1, 0, 3, 2, 1};
1064ef8b56bfSJed Brown   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};
10659371c9d4SSatish Balay   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};
1066ef8b56bfSJed Brown   static const PetscInt tripInv[12]  = {-5, -6, -4, -3, -2, -1, 0, 2, 1, 3, 4, 5};
1067ef8b56bfSJed Brown   static const PetscInt ttriInv[12]  = {-6, -5, -4, -3, -2, -1, 0, 2, 1, 3, 5, 4};
1068ef8b56bfSJed Brown   static const PetscInt tquadInv[16] = {-8, -7, -6, -5, -4, -3, -2, -1, 0, 3, 2, 1, 4, 7, 6, 5};
1069ef8b56bfSJed Brown   static const PetscInt pyrInv[8]    = {-4, -3, -2, -1, 0, 3, 2, 1};
1070b5a892a1SMatthew G. Knepley   switch (ct) {
1071d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
1072d71ae5a4SJacob Faibussowitsch     return 0;
1073b5a892a1SMatthew G. Knepley   case DM_POLYTOPE_SEGMENT:
1074d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
1075d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, o2);
1076d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
1077d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, triInv[o2 + 3]);
1078d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
1079d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, quadInv[o2 + 4]);
1080d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
1081d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, o2);
1082d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
1083d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, tetInv[o2 + 12]);
1084d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
1085d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, hexInv[o2 + 24]);
1086d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
1087d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, tripInv[o2 + 6]);
1088d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
1089d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, ttriInv[o2 + 6]);
1090d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
1091d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, tquadInv[o2 + 8]);
1092d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
1093d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, pyrInv[o2 + 4]);
1094d71ae5a4SJacob Faibussowitsch   default:
1095d71ae5a4SJacob Faibussowitsch     return 0;
1096b5a892a1SMatthew G. Knepley   }
1097b5a892a1SMatthew G. Knepley }
1098b5a892a1SMatthew G. Knepley 
1099b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeMatchOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *);
1100b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeMatchVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *);
1101b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeGetOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *);
1102b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeGetVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *);
1103012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeInCellTest(DMPolytopeType, const PetscReal[], PetscBool *);
1104