xref: /petsc/include/petscdm.h (revision 974ca4ec8a1e8cc5f700b20b09fb20d3d9f6c85d)
1e1589f56SBarry Smith /*
2e1589f56SBarry Smith       Objects to manage the interactions between the mesh data structures and the algebraic objects
3e1589f56SBarry Smith */
46524c165SJacob Faibussowitsch #ifndef PETSCDM_H
526bd1501SBarry Smith #define PETSCDM_H
62c8e378dSBarry Smith #include <petscmat.h>
71e25c274SJed Brown #include <petscdmtypes.h>
8decb47aaSMatthew G. Knepley #include <petscfetypes.h>
92764a2aaSMatthew G. Knepley #include <petscdstypes.h>
10c58f1c22SToby Isaac #include <petscdmlabel.h>
1107218a29SMatthew G. Knepley #include <petscdt.h>
12e1589f56SBarry Smith 
13ac09b921SBarry Smith /* SUBMANSEC = DM */
14ac09b921SBarry Smith 
15607a6623SBarry Smith PETSC_EXTERN PetscErrorCode DMInitializePackage(void);
16e1589f56SBarry Smith 
17014dd563SJed Brown PETSC_EXTERN PetscClassId DM_CLASSID;
18e1589f56SBarry Smith 
19528c63e0SDave May #define DMLOCATEPOINT_POINT_NOT_FOUND -367
20528c63e0SDave May 
2176bdecfbSBarry Smith /*J
2287497f52SBarry Smith     DMType - String with the name of a PETSc `DM`
23e1589f56SBarry Smith 
24e1589f56SBarry Smith    Level: beginner
25e1589f56SBarry Smith 
261cc06b55SBarry Smith .seealso: [](ch_dmbase), `DMSetType()`, `DMCreate()`, `DM`
2776bdecfbSBarry Smith J*/
2819fd82e9SBarry Smith typedef const char *DMType;
29e1589f56SBarry Smith #define DMDA        "da"
30e1589f56SBarry Smith #define DMCOMPOSITE "composite"
31e1589f56SBarry Smith #define DMSLICED    "sliced"
32fe1899a2SJed Brown #define DMSHELL     "shell"
33ab7f58a0SBarry Smith #define DMPLEX      "plex"
348ac4e037SJed Brown #define DMREDUNDANT "redundant"
353a19ef87SMatthew G Knepley #define DMPATCH     "patch"
361d72bce8STim Tautges #define DMMOAB      "moab"
375f2c45f1SShri Abhyankar #define DMNETWORK   "network"
38ef51cf95SToby Isaac #define DMFOREST    "forest"
39db4d5e8cSToby Isaac #define DMP4EST     "p4est"
40db4d5e8cSToby Isaac #define DMP8EST     "p8est"
412fd35b1fSDave May #define DMSWARM     "swarm"
42d852a638SPatrick Sanan #define DMPRODUCT   "product"
43a3101111SPatrick Sanan #define DMSTAG      "stag"
44e1589f56SBarry Smith 
45bff4a2f0SMatthew G. Knepley PETSC_EXTERN const char *const       DMBoundaryTypes[];
4640967b3bSMatthew G. Knepley PETSC_EXTERN const char *const       DMBoundaryConditionTypes[];
47863027abSJed Brown PETSC_EXTERN const char *const       DMBlockingTypes[];
48140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList       DMList;
49c0517cd5SMatthew G. Knepley PETSC_EXTERN DMGeneratorFunctionList DMGenerateList;
50014dd563SJed Brown PETSC_EXTERN PetscErrorCode          DMCreate(MPI_Comm, DM *);
5138221697SMatthew G. Knepley PETSC_EXTERN PetscErrorCode          DMClone(DM, DM *);
5219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode          DMSetType(DM, DMType);
5319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode          DMGetType(DM, DMType *);
54bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode          DMRegister(const char[], PetscErrorCode (*)(DM));
55014dd563SJed Brown PETSC_EXTERN PetscErrorCode          DMRegisterDestroy(void);
56e1589f56SBarry Smith 
57014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMView(DM, PetscViewer);
58014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLoad(DM, PetscViewer);
59014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDestroy(DM *);
60014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM, Vec *);
61014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM, Vec *);
62014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM, Vec *);
63014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM, Vec *);
64014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM, Vec *);
65014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM, Vec *);
66014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
6750eeb1caSToby Isaac PETSC_EXTERN PetscErrorCode DMClearLocalVectors(DM);
68*974ca4ecSStefano Zampini PETSC_EXTERN PetscErrorCode DMClearNamedGlobalVectors(DM);
69*974ca4ecSStefano Zampini PETSC_EXTERN PetscErrorCode DMClearNamedLocalVectors(DM);
70e77ac854SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasNamedGlobalVector(DM, const char *, PetscBool *);
71014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM, const char *, Vec *);
72014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM, const char *, Vec *);
73e77ac854SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasNamedLocalVector(DM, const char *, PetscBool *);
742348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM, const char *, Vec *);
752348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM, const char *, Vec *);
76014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM, ISLocalToGlobalMapping *);
77014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM, PetscInt *, char ***, IS **);
78014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM, PetscInt *);
79b412c318SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateColoring(DM, ISColoringType, ISColoring *);
80b412c318SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM, Mat *);
81aa0f6e3cSJed Brown PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateSkip(DM, PetscBool);
82014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM, PetscBool);
83b06ff27eSHong Zhang PETSC_EXTERN PetscErrorCode DMSetMatrixStructureOnly(DM, PetscBool);
84863027abSJed Brown PETSC_EXTERN PetscErrorCode DMSetBlockingType(DM, DMBlockingType);
85863027abSJed Brown PETSC_EXTERN PetscErrorCode DMGetBlockingType(DM, DMBlockingType *);
86014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM, DM, Mat *, Vec *);
873ad4599aSBarry Smith PETSC_EXTERN PetscErrorCode DMCreateRestriction(DM, DM, Mat *);
886dbf9973SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMCreateInjection(DM, DM, Mat *);
89bd041c0cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateMassMatrix(DM, DM, Mat *);
90b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateMassMatrixLumped(DM, Vec *);
9169291d52SBarry Smith PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM, PetscInt, MPI_Datatype, void *);
9269291d52SBarry Smith PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM, PetscInt, MPI_Datatype, void *);
93014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefine(DM, MPI_Comm, DM *);
94014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsen(DM, MPI_Comm, DM *);
95a8fb8f29SToby Isaac PETSC_EXTERN PetscErrorCode DMGetCoarseDM(DM, DM *);
96a8fb8f29SToby Isaac PETSC_EXTERN PetscErrorCode DMSetCoarseDM(DM, DM);
9788bdff64SToby Isaac PETSC_EXTERN PetscErrorCode DMGetFineDM(DM, DM *);
9888bdff64SToby Isaac PETSC_EXTERN PetscErrorCode DMSetFineDM(DM, DM);
99014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM, PetscInt, DM[]);
100014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM, PetscInt, DM[]);
101014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, Vec, Mat, DM, void *), void *);
102dc822a44SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, Vec, Mat, DM, void *), void *);
103014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, DM, void *), void *);
1043d8e3701SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, DM, void *), void *);
105014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestrict(DM, Mat, Vec, Mat, DM);
106014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMInterpolate(DM, Mat, DM);
1071f3379b2SToby Isaac PETSC_EXTERN PetscErrorCode DMInterpolateSolution(DM, DM, Mat, Vec, Vec);
108d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMExtrude(DM, PetscInt, DM *);
109014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
110fe2efc57SMark PETSC_EXTERN PetscErrorCode DMViewFromOptions(DM, PetscObject, const char[]);
111ca266f36SBarry Smith 
112c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerate(DM, const char[], PetscBool, DM *);
1139fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMGenerateRegister(const char[], PetscErrorCode (*)(DM, PetscBool, DM *), PetscErrorCode (*)(DM, PetscReal *, DM *), PetscErrorCode (*)(DM, Vec, DMLabel, DMLabel, DM *), PetscInt);
114c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerateRegisterAll(void);
115c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerateRegisterDestroy(void);
116a1b0c543SToby Isaac PETSC_EXTERN PetscErrorCode DMAdaptLabel(DM, DMLabel, DM *);
1179fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMAdaptMetric(DM, Vec, DMLabel, DMLabel, DM *);
118df0b854cSToby Isaac 
119014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetUp(DM);
120014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM, DM, Mat, Vec *);
12197779f9aSLisandro Dalcin PETSC_EXTERN                PETSC_DEPRECATED_FUNCTION("Use DMDACreateAggregates() or DMCreateRestriction() (since version 3.12)") PetscErrorCode DMCreateAggregates(DM, DM, Mat *);
122baf369e7SPeter Brune PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM, PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), void *);
123d4d07f1eSToby Isaac PETSC_EXTERN PetscErrorCode DMLocalToGlobalHookAdd(DM, PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), void *);
12401729b5cSPatrick Sanan PETSC_EXTERN PetscErrorCode DMGlobalToLocal(DM, Vec, InsertMode, Vec);
125014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM, Vec, InsertMode, Vec);
126014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM, Vec, InsertMode, Vec);
12701729b5cSPatrick Sanan PETSC_EXTERN PetscErrorCode DMLocalToGlobal(DM, Vec, InsertMode, Vec);
128014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM, Vec, InsertMode, Vec);
129014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM, Vec, InsertMode, Vec);
130d78e899eSRichard Tran Mills PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM, Vec, InsertMode, Vec);
131d78e899eSRichard Tran Mills PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM, Vec, InsertMode, Vec);
13219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMConvert(DM, DMType, DM *);
133e1589f56SBarry Smith 
134c73cfb54SMatthew G. Knepley /* Topology support */
135c73cfb54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDimension(DM, PetscInt *);
136c73cfb54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetDimension(DM, PetscInt);
137793f3fe5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDimPoints(DM, PetscInt, PetscInt *, PetscInt *);
1388e4ac7eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetUseNatural(DM, PetscBool *);
1398e4ac7eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetUseNatural(DM, PetscBool);
1406858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNeighbors(DM, PetscInt *, const PetscMPIInt **);
14146e270d4SMatthew G. Knepley 
14246e270d4SMatthew G. Knepley /* Coordinate support */
1436636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM, DM *);
1441cfe2091SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM, DM);
1456858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinateDM(DM, DM *);
1466858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinateDM(DM, DM);
14746e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateDim(DM, PetscInt *);
14846e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateDim(DM, PetscInt);
149e8abe2deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM, PetscSection *);
15046e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM, PetscInt, PetscSection);
1516858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinateSection(DM, PetscSection *);
1526858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinateSection(DM, PetscInt, PetscSection);
1536636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM, Vec *);
1546636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM, Vec);
1556858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinates(DM, Vec *);
1566858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinates(DM, Vec);
15781e9a530SVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalSetUp(DM);
1586858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM, Vec *);
15981e9a530SVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalNoncollective(DM, Vec *);
1602db98f8dSVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalTuple(DM, IS, PetscSection *, Vec *);
1616636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM, Vec);
1626858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM);
1636858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocal(DM, Vec *);
1646858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM, Vec *);
1656858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinatesLocal(DM, Vec);
1666858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateField(DM, DMField *);
1676858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateField(DM, DMField);
1686858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetLocalBoundingBox(DM, PetscReal[], PetscReal[]);
1696858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetBoundingBox(DM, PetscReal[], PetscReal[]);
1706858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectCoordinates(DM, PetscFE);
17162a38674SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocatePoints(DM, Vec, DMPointLocationType, PetscSF *);
1726858538eSMatthew G. Knepley 
1736858538eSMatthew G. Knepley /* Periodicity support */
1744fb89dddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM, const PetscReal *[], const PetscReal *[], const PetscReal *[]);
1754fb89dddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM, const PetscReal[], const PetscReal[], const PetscReal[]);
176e907e85cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate(DM, const PetscScalar[], PetscBool, PetscScalar[]);
1772e17dfb7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocalizeCoordinates(DM);
17836447a5eSToby Isaac PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalized(DM, PetscBool *);
1798f700142SStefano Zampini PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalizedLocal(DM, PetscBool *);
1806636e97aSMatthew G Knepley 
1815dbd56e3SPeter Brune /* block hook interface */
182be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, VecScatter, VecScatter, DM, void *), void *);
183b3a6b972SJed Brown PETSC_EXTERN PetscErrorCode DMSubDomainHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, VecScatter, VecScatter, DM, void *), void *);
184be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM, VecScatter, VecScatter, DM);
1855dbd56e3SPeter Brune 
186014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM, const char[]);
18731697293SDave May PETSC_EXTERN PetscErrorCode DMAppendOptionsPrefix(DM, const char[]);
18831697293SDave May PETSC_EXTERN PetscErrorCode DMGetOptionsPrefix(DM, const char *[]);
18919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetVecType(DM, VecType);
190c0dedaeaSBarry Smith PETSC_EXTERN PetscErrorCode DMGetVecType(DM, VecType *);
19119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetMatType(DM, MatType);
192c0dedaeaSBarry Smith PETSC_EXTERN PetscErrorCode DMGetMatType(DM, MatType *);
1938f1509bcSBarry Smith PETSC_EXTERN PetscErrorCode DMSetISColoringType(DM, ISColoringType);
1948f1509bcSBarry Smith PETSC_EXTERN PetscErrorCode DMGetISColoringType(DM, ISColoringType *);
195014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM, void *);
196014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM, PetscErrorCode (*)(void **));
197014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM, void *);
198014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM, PetscErrorCode (*)(DM, Vec, Vec));
199014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM, PetscBool *);
200b0ae01b7SPeter Brune PETSC_EXTERN PetscErrorCode DMHasColoring(DM, PetscBool *);
2013ad4599aSBarry Smith PETSC_EXTERN PetscErrorCode DMHasCreateRestriction(DM, PetscBool *);
202a7058e45SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMHasCreateInjection(DM, PetscBool *);
203014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM, Vec, Vec);
20493d92d96SBarry Smith 
20537bc7515SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, const PetscInt[], IS *, DM *);
2062adcc780SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSuperDM(DM[], PetscInt, IS **, DM *);
207792b654fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSectionSubDM(DM, PetscInt, const PetscInt[], IS *, DM *);
208792b654fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSectionSuperDM(DM[], PetscInt, IS **, DM *);
20916621825SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM, PetscInt *, char ***, IS **, DM **);
2108d4ac253SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM, PetscInt *, char ***, IS **, IS **, DM **);
211e30e807fSPeter Brune PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM, PetscInt, DM *, VecScatter **, VecScatter **, VecScatter **);
212e7c4fc90SDmitry Karpeev 
213014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM, PetscInt *);
214fef3a512SBarry Smith PETSC_EXTERN PetscErrorCode DMSetRefineLevel(DM, PetscInt);
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM, PetscInt *);
2169a64c4a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoarsenLevel(DM, PetscInt);
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);
218e1589f56SBarry Smith 
2195f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM *);
2205f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
221c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM *);
222c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
223531c7667SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringUseDM(Mat, MatFDColoring);
2245f1ad066SMatthew G Knepley 
225e1589f56SBarry Smith typedef struct NLF_DAAD *NLF;
226e1589f56SBarry Smith 
227bc2bf880SBarry Smith #define DM_FILE_CLASSID 1211221
2287da65231SMatthew G Knepley 
2297da65231SMatthew G Knepley /* FEM support */
230014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char[], PetscInt, const PetscScalar[]);
231014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char[], PetscInt, PetscInt, const PetscScalar[]);
2326113b454SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char[], PetscReal, Vec);
2337da65231SMatthew G Knepley 
2348cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *));
2358cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *));
2368cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *));
2378cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *));
238f9d4088aSMatthew G. Knepley 
239061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMGetSection(DM, PetscSection *); /* Use DMGetLocalSection() in new code (since v3.12) */
240061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMSetSection(DM, PetscSection);   /* Use DMSetLocalSection() in new code (since v3.12) */
241061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalSection(DM, PetscSection *);
242061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMSetLocalSection(DM, PetscSection);
243e87a4003SBarry Smith PETSC_EXTERN PetscErrorCode DMGetGlobalSection(DM, PetscSection *);
244e87a4003SBarry Smith PETSC_EXTERN PetscErrorCode DMSetGlobalSection(DM, PetscSection);
245d71ae5a4SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Use DMGetSection() (since v3.9)") PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *s)
246d71ae5a4SJacob Faibussowitsch {
2479371c9d4SSatish Balay   return DMGetSection(dm, s);
2489371c9d4SSatish Balay }
249d71ae5a4SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Use DMSetSection() (since v3.9)") PetscErrorCode DMSetDefaultSection(DM dm, PetscSection s)
250d71ae5a4SJacob Faibussowitsch {
2519371c9d4SSatish Balay   return DMSetSection(dm, s);
2529371c9d4SSatish Balay }
253d71ae5a4SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Use DMGetGlobalSection() (since v3.9)") PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *s)
254d71ae5a4SJacob Faibussowitsch {
2559371c9d4SSatish Balay   return DMGetGlobalSection(dm, s);
2569371c9d4SSatish Balay }
257d71ae5a4SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Use DMSetGlobalSection() (since v3.9)") PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection s)
258d71ae5a4SJacob Faibussowitsch {
2599371c9d4SSatish Balay   return DMSetGlobalSection(dm, s);
2609371c9d4SSatish Balay }
261e87a4003SBarry Smith 
2621bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMGetSectionSF(DM, PetscSF *);
2631bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMSetSectionSF(DM, PetscSF);
2641bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateSectionSF(DM, PetscSection, PetscSection);
265d71ae5a4SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Use DMGetSectionSF() (since v3.12)") PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *s)
266d71ae5a4SJacob Faibussowitsch {
2679371c9d4SSatish Balay   return DMGetSectionSF(dm, s);
2689371c9d4SSatish Balay }
269d71ae5a4SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Use DMSetSectionSF() (since v3.12)") PetscErrorCode DMSetDefaultSF(DM dm, PetscSF s)
270d71ae5a4SJacob Faibussowitsch {
2719371c9d4SSatish Balay   return DMSetSectionSF(dm, s);
2729371c9d4SSatish Balay }
273d71ae5a4SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Use DMCreateSectionSF() (since v3.12)") PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection l, PetscSection g)
274d71ae5a4SJacob Faibussowitsch {
2759371c9d4SSatish Balay   return DMCreateSectionSF(dm, l, g);
2769371c9d4SSatish Balay }
2771bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
2781bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);
2794f37162bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNaturalSF(DM, PetscSF *);
2804f37162bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNaturalSF(DM, PetscSF);
2811bb6d2a8SBarry Smith 
28279769bd5SJed Brown PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *, Vec *);
28379769bd5SJed Brown PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat, Vec);
284e87a4003SBarry Smith 
28514f150ffSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *);
286cdb7a50dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *);
287cdb7a50dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal);
288cdb7a50dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char *, PetscInt, PetscReal *);
28914f150ffSMatthew G. Knepley 
290af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
291af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
29244a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, DMLabel *, PetscObject *);
29344a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, DMLabel, PetscObject);
29444a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAddField(DM, DMLabel, PetscObject);
295e0b68406SMatthew Knepley PETSC_EXTERN PetscErrorCode DMSetFieldAvoidTensor(DM, PetscInt, PetscBool);
296e0b68406SMatthew Knepley PETSC_EXTERN PetscErrorCode DMGetFieldAvoidTensor(DM, PetscInt, PetscBool *);
29744a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMClearFields(DM);
298e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyFields(DM, DM);
29934aa8a36SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAdjacency(DM, PetscInt, PetscBool *, PetscBool *);
30034aa8a36SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetAdjacency(DM, PetscInt, PetscBool, PetscBool);
301b0441da4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetBasicAdjacency(DM, PetscBool *, PetscBool *);
302b0441da4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetBasicAdjacency(DM, PetscBool, PetscBool);
303e5e52638SMatthew G. Knepley 
304e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNumDS(DM, PetscInt *);
305e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *);
30607218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellDS(DM, PetscInt, PetscDS *, PetscDS *);
30707218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetRegionDS(DM, DMLabel, IS *, PetscDS *, PetscDS *);
30807218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetRegionDS(DM, DMLabel, IS, PetscDS, PetscDS);
30907218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetRegionNumDS(DM, PetscInt, DMLabel *, IS *, PetscDS *, PetscDS *);
31007218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetRegionNumDS(DM, PetscInt, DMLabel, IS, PetscDS, PetscDS);
3111d3af9e0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMFindRegionNum(DM, PetscDS, PetscInt *);
3122df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateFEDefault(DM, PetscInt, const char[], PetscInt, PetscFE *);
313e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateDS(DM);
314e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMClearDS(DM);
315e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyDS(DM, DM);
316e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyDisc(DM, DM);
317f2cacb80SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMComputeExactSolution(DM, PetscReal, Vec, Vec);
3189a2a23afSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNumAuxiliaryVec(DM, PetscInt *);
319ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAuxiliaryVec(DM, DMLabel, PetscInt, PetscInt, Vec *);
320ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetAuxiliaryVec(DM, DMLabel, PetscInt, PetscInt, Vec);
321ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAuxiliaryLabels(DM, DMLabel[], PetscInt[], PetscInt[]);
3229a2a23afSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyAuxiliaryVec(DM, DM);
323af122d2aSMatthew G Knepley 
3244267b1a3SMatthew G. Knepley /*MC
3254267b1a3SMatthew G. Knepley   DMInterpolationInfo - Structure for holding information about interpolation on a mesh
3264267b1a3SMatthew G. Knepley 
3274267b1a3SMatthew G. Knepley   Synopsis:
3284267b1a3SMatthew G. Knepley     comm   - The communicator
3294267b1a3SMatthew G. Knepley     dim    - The spatial dimension of points
3304267b1a3SMatthew G. Knepley     nInput - The number of input points
3314267b1a3SMatthew G. Knepley     points - The input point coordinates
3324267b1a3SMatthew G. Knepley     cells  - The cell containing each point
3334267b1a3SMatthew G. Knepley     n      - The number of local points
3344267b1a3SMatthew G. Knepley     coords - The point coordinates
3354267b1a3SMatthew G. Knepley     dof    - The number of components to interpolate
3364267b1a3SMatthew G. Knepley 
33716a05f60SBarry Smith   Level: intermediate
33816a05f60SBarry Smith 
33916a05f60SBarry Smith .seealso: `DM`, `DMInterpolationCreate()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
3404267b1a3SMatthew G. Knepley M*/
341e87bb0d3SMatthew G Knepley struct _DMInterpolationInfo {
342e87bb0d3SMatthew G Knepley   MPI_Comm   comm;
343e87bb0d3SMatthew G Knepley   PetscInt   dim;    /*1 The spatial dimension of points */
344e87bb0d3SMatthew G Knepley   PetscInt   nInput; /* The number of input points */
345e87bb0d3SMatthew G Knepley   PetscReal *points; /* The input point coordinates */
346e87bb0d3SMatthew G Knepley   PetscInt  *cells;  /* The cell containing each point */
347e87bb0d3SMatthew G Knepley   PetscInt   n;      /* The number of local points */
348e87bb0d3SMatthew G Knepley   Vec        coords; /* The point coordinates */
349e87bb0d3SMatthew G Knepley   PetscInt   dof;    /* The number of components to interpolate */
350e87bb0d3SMatthew G Knepley };
351e87bb0d3SMatthew G Knepley typedef struct _DMInterpolationInfo *DMInterpolationInfo;
352e87bb0d3SMatthew G Knepley 
35394b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
35494b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
35594b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
35694b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
35794b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
35894b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
35952aa1562SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool, PetscBool);
36094b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
36194b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
36294b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
36394b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
36494b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);
365c58f1c22SToby Isaac 
366c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMCreateLabel(DM, const char[]);
367c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelValue(DM, const char[], PetscInt, PetscInt *);
368c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMSetLabelValue(DM, const char[], PetscInt, PetscInt);
369c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMClearLabelValue(DM, const char[], PetscInt, PetscInt);
370c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelSize(DM, const char[], PetscInt *);
371c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelIdIS(DM, const char[], IS *);
372c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetStratumSize(DM, const char[], PetscInt, PetscInt *);
373c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetStratumIS(DM, const char[], PetscInt, IS *);
3744de306b1SToby Isaac PETSC_EXTERN PetscErrorCode DMSetStratumIS(DM, const char[], PetscInt, IS);
375c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMClearLabelStratum(DM, const char[], PetscInt);
376c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelOutput(DM, const char[], PetscBool *);
377c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMSetLabelOutput(DM, const char[], PetscBool);
3785f15299fSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMGetFirstLabeledPoint(DM, DM, DMLabel, PetscInt, const PetscInt *, PetscInt, PetscInt *, PetscDS *);
379c58f1c22SToby Isaac 
3802cbb9b06SVaclav Hapla /*E
38187497f52SBarry Smith    DMCopyLabelsMode - Determines how `DMCopyLabels()` behaves when there is a `DMLabel` in the source and destination DMs with the same name
3822cbb9b06SVaclav Hapla 
38316a05f60SBarry Smith    Values:
38416a05f60SBarry Smith +  `DM_COPY_LABELS_REPLACE`  - replace label in destination by label from source
38516a05f60SBarry Smith .  `DM_COPY_LABELS_KEEP`     - keep destination label
38616a05f60SBarry Smith -  `DM_COPY_LABELS_FAIL`     - throw error
38716a05f60SBarry Smith 
3882cbb9b06SVaclav Hapla    Level: advanced
3892cbb9b06SVaclav Hapla 
39016a05f60SBarry Smith .seealso: `DMLabel`, `DM`, `DMCompareLabels()`, `DMRemoveLabel()`
3912cbb9b06SVaclav Hapla E*/
3929371c9d4SSatish Balay typedef enum {
3939371c9d4SSatish Balay   DM_COPY_LABELS_REPLACE,
3949371c9d4SSatish Balay   DM_COPY_LABELS_KEEP,
3959371c9d4SSatish Balay   DM_COPY_LABELS_FAIL
3969371c9d4SSatish Balay } DMCopyLabelsMode;
3972cbb9b06SVaclav Hapla PETSC_EXTERN const char *const DMCopyLabelsModes[];
3982cbb9b06SVaclav Hapla 
399c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetNumLabels(DM, PetscInt *);
400c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelName(DM, PetscInt, const char **);
401c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMHasLabel(DM, const char[], PetscBool *);
402c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabel(DM, const char *, DMLabel *);
4034a7ee7d0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetLabel(DM, DMLabel);
404c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelByNum(DM, PetscInt, DMLabel *);
405c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMAddLabel(DM, DMLabel);
406c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMRemoveLabel(DM, const char[], DMLabel *);
407306894acSVaclav Hapla PETSC_EXTERN PetscErrorCode DMRemoveLabelBySelf(DM, DMLabel *, PetscBool);
4082cbb9b06SVaclav Hapla PETSC_EXTERN PetscErrorCode DMCopyLabels(DM, DM, PetscCopyMode, PetscBool, DMCopyLabelsMode emode);
409609dae6eSVaclav Hapla PETSC_EXTERN PetscErrorCode DMCompareLabels(DM, DM, PetscBool *, char **);
410c58f1c22SToby Isaac 
41145480ffeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAddBoundary(DM, DMBoundaryConditionType, const char[], DMLabel, PetscInt, const PetscInt[], PetscInt, PetscInt, const PetscInt[], void (*)(void), void (*)(void), void *, PetscInt *);
412a6ba4734SToby Isaac PETSC_EXTERN PetscErrorCode DMIsBoundaryPoint(DM, PetscInt, PetscBool *);
4134d6f44ffSToby Isaac 
4140709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMProjectFunction(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
4150709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMProjectFunctionLocal(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
4162c53366bSMatthew 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);
4171c531cf8SMatthew 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);
418191494d9SMatthew 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);
419d29d7c6eSMatthew 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);
4201c531cf8SMatthew 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);
421ece3a9fcSMatthew 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);
4220709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2Diff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
423b698f381SToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2GradientDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal[], PetscReal *);
4241189c1efSToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2FieldDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
4252e4af2aeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMComputeError(DM, Vec, PetscReal[], Vec *);
426ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasBasisTransform(DM, PetscBool *);
427ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyTransform(DM, DM);
4288320bc6fSPatrick Sanan 
4298320bc6fSPatrick Sanan PETSC_EXTERN PetscErrorCode DMGetCompatibility(DM, DM, PetscBool *, PetscBool *);
430c0f0dcc3SMatthew G. Knepley 
431c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitorSet(DM, PetscErrorCode (*)(DM, void *), void *, PetscErrorCode (*)(void **));
432c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitorCancel(DM);
433c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitorSetFromOptions(DM, const char[], const char[], const char[], PetscErrorCode (*)(DM, void *), PetscErrorCode (*)(DM, PetscViewerAndFormat *), PetscBool *);
434c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitor(DM);
435c0f0dcc3SMatthew G. Knepley 
436d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetDim(DMPolytopeType ct)
437d71ae5a4SJacob Faibussowitsch {
438412e9a14SMatthew G. Knepley   switch (ct) {
439d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
440d71ae5a4SJacob Faibussowitsch     return 0;
441412e9a14SMatthew G. Knepley   case DM_POLYTOPE_SEGMENT:
442d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
443d71ae5a4SJacob Faibussowitsch     return 1;
444412e9a14SMatthew G. Knepley   case DM_POLYTOPE_TRIANGLE:
445412e9a14SMatthew G. Knepley   case DM_POLYTOPE_QUADRILATERAL:
446d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
447d71ae5a4SJacob Faibussowitsch     return 2;
448412e9a14SMatthew G. Knepley   case DM_POLYTOPE_TETRAHEDRON:
449412e9a14SMatthew G. Knepley   case DM_POLYTOPE_HEXAHEDRON:
450412e9a14SMatthew G. Knepley   case DM_POLYTOPE_TRI_PRISM:
451412e9a14SMatthew G. Knepley   case DM_POLYTOPE_TRI_PRISM_TENSOR:
452412e9a14SMatthew G. Knepley   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
453d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
454d71ae5a4SJacob Faibussowitsch     return 3;
455d71ae5a4SJacob Faibussowitsch   default:
456d71ae5a4SJacob Faibussowitsch     return -1;
457412e9a14SMatthew G. Knepley   }
458412e9a14SMatthew G. Knepley }
459412e9a14SMatthew G. Knepley 
460d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetConeSize(DMPolytopeType ct)
461d71ae5a4SJacob Faibussowitsch {
462412e9a14SMatthew G. Knepley   switch (ct) {
463d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
464d71ae5a4SJacob Faibussowitsch     return 0;
465d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
466d71ae5a4SJacob Faibussowitsch     return 2;
467d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
468d71ae5a4SJacob Faibussowitsch     return 2;
469d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
470d71ae5a4SJacob Faibussowitsch     return 3;
471d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
472d71ae5a4SJacob Faibussowitsch     return 4;
473d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
474d71ae5a4SJacob Faibussowitsch     return 4;
475d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
476d71ae5a4SJacob Faibussowitsch     return 4;
477d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
478d71ae5a4SJacob Faibussowitsch     return 6;
479d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
480d71ae5a4SJacob Faibussowitsch     return 5;
481d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
482d71ae5a4SJacob Faibussowitsch     return 5;
483d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
484d71ae5a4SJacob Faibussowitsch     return 6;
485d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
486d71ae5a4SJacob Faibussowitsch     return 5;
487d71ae5a4SJacob Faibussowitsch   default:
488d71ae5a4SJacob Faibussowitsch     return -1;
489412e9a14SMatthew G. Knepley   }
490412e9a14SMatthew G. Knepley }
491412e9a14SMatthew G. Knepley 
492d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetNumVertices(DMPolytopeType ct)
493d71ae5a4SJacob Faibussowitsch {
494412e9a14SMatthew G. Knepley   switch (ct) {
495d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
496d71ae5a4SJacob Faibussowitsch     return 1;
497d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
498d71ae5a4SJacob Faibussowitsch     return 2;
499d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
500d71ae5a4SJacob Faibussowitsch     return 2;
501d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
502d71ae5a4SJacob Faibussowitsch     return 3;
503d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
504d71ae5a4SJacob Faibussowitsch     return 4;
505d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
506d71ae5a4SJacob Faibussowitsch     return 4;
507d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
508d71ae5a4SJacob Faibussowitsch     return 4;
509d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
510d71ae5a4SJacob Faibussowitsch     return 8;
511d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
512d71ae5a4SJacob Faibussowitsch     return 6;
513d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
514d71ae5a4SJacob Faibussowitsch     return 6;
515d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
516d71ae5a4SJacob Faibussowitsch     return 8;
517d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
518d71ae5a4SJacob Faibussowitsch     return 5;
519d71ae5a4SJacob Faibussowitsch   default:
520d71ae5a4SJacob Faibussowitsch     return -1;
521412e9a14SMatthew G. Knepley   }
522412e9a14SMatthew G. Knepley }
523412e9a14SMatthew G. Knepley 
524d71ae5a4SJacob Faibussowitsch static inline DMPolytopeType DMPolytopeTypeSimpleShape(PetscInt dim, PetscBool simplex)
525d71ae5a4SJacob Faibussowitsch {
5269371c9d4SSatish 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)));
5279318fe57SMatthew G. Knepley }
5289318fe57SMatthew G. Knepley 
529d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetNumArrangments(DMPolytopeType ct)
530d71ae5a4SJacob Faibussowitsch {
531b5a892a1SMatthew G. Knepley   switch (ct) {
532d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
533d71ae5a4SJacob Faibussowitsch     return 1;
534d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
535d71ae5a4SJacob Faibussowitsch     return 2;
536d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
537d71ae5a4SJacob Faibussowitsch     return 2;
538d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
539d71ae5a4SJacob Faibussowitsch     return 6;
540d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
541d71ae5a4SJacob Faibussowitsch     return 8;
542d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
543d71ae5a4SJacob Faibussowitsch     return 4;
544d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
545d71ae5a4SJacob Faibussowitsch     return 24;
546d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
547d71ae5a4SJacob Faibussowitsch     return 48;
548d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
549d71ae5a4SJacob Faibussowitsch     return 12;
550d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
551d71ae5a4SJacob Faibussowitsch     return 12;
552d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
553d71ae5a4SJacob Faibussowitsch     return 16;
554d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
555d71ae5a4SJacob Faibussowitsch     return 8;
556d71ae5a4SJacob Faibussowitsch   default:
557d71ae5a4SJacob Faibussowitsch     return -1;
558b5a892a1SMatthew G. Knepley   }
559b5a892a1SMatthew G. Knepley }
560b5a892a1SMatthew G. Knepley 
561b5a892a1SMatthew G. Knepley /* An arrangement is a face order combined with an orientation for each face */
562d71ae5a4SJacob Faibussowitsch static inline const PetscInt *DMPolytopeTypeGetArrangment(DMPolytopeType ct, PetscInt o)
563d71ae5a4SJacob Faibussowitsch {
564ef8b56bfSJed Brown   static const PetscInt pntArr[1 * 2] = {0, 0};
565b5a892a1SMatthew G. Knepley   /* a: swap */
5669371c9d4SSatish Balay   static const PetscInt segArr[2 * 2 * 2] = {1, 0, 0, 0, /* -1: a */
5679371c9d4SSatish Balay                                              0, 0, 1, 0,
5689371c9d4SSatish Balay                                              /*  0: e */};
569b5a892a1SMatthew G. Knepley   /* a: swap first two
570b5a892a1SMatthew G. Knepley      b: swap last two */
5719371c9d4SSatish Balay   static const PetscInt triArr[6 * 3 * 2] = {0, -1, 2, -1, 1, -1, /* -3: b */
572b5a892a1SMatthew G. Knepley                                              2, -1, 1, -1, 0, -1, /* -2: aba */
573b5a892a1SMatthew G. Knepley                                              1, -1, 0, -1, 2, -1, /* -1: a */
574b5a892a1SMatthew G. Knepley                                              0, 0,  1, 0,  2, 0,  /*  0: identity */
575b5a892a1SMatthew G. Knepley                                              1, 0,  2, 0,  0, 0,  /*  1: ba */
5769371c9d4SSatish Balay                                              2, 0,  0, 0,  1, 0,
5779371c9d4SSatish Balay                                              /*  2: ab */};
578b5a892a1SMatthew G. Knepley   /* a: forward cyclic permutation
579b5a892a1SMatthew G. Knepley      b: swap first and last pairs */
5809371c9d4SSatish Balay   static const PetscInt quadArr[8 * 4 * 2] = {1, -1, 0, -1, 3, -1, 2, -1, /* -4: b */
581b5a892a1SMatthew G. Knepley                                               0, -1, 3, -1, 2, -1, 1, -1, /* -3: b a^3 = a b */
582b5a892a1SMatthew G. Knepley                                               3, -1, 2, -1, 1, -1, 0, -1, /* -2: b a^2 = a^2 b */
583b5a892a1SMatthew G. Knepley                                               2, -1, 1, -1, 0, -1, 3, -1, /* -1: b a   = a^3 b */
584b5a892a1SMatthew G. Knepley                                               0, 0,  1, 0,  2, 0,  3, 0,  /*  0: identity */
585b5a892a1SMatthew G. Knepley                                               1, 0,  2, 0,  3, 0,  0, 0,  /*  1: a */
586b5a892a1SMatthew G. Knepley                                               2, 0,  3, 0,  0, 0,  1, 0,  /*  2: a^2 */
5879371c9d4SSatish Balay                                               3, 0,  0, 0,  1, 0,  2, 0,
5889371c9d4SSatish Balay                                               /*  3: a^3 */};
589b5a892a1SMatthew G. Knepley   /* r: rotate 180
590b5a892a1SMatthew G. Knepley      b: swap top and bottom segments */
5919371c9d4SSatish Balay   static const PetscInt tsegArr[4 * 4 * 2] = {1, -1, 0, -1, 3, -1, 2, -1, /* -2: r b */
592b5a892a1SMatthew G. Knepley                                               0, -1, 1, -1, 3, 0,  2, 0,  /* -1: r */
593b5a892a1SMatthew G. Knepley                                               0, 0,  1, 0,  2, 0,  3, 0,  /*  0: identity */
5949371c9d4SSatish Balay                                               1, 0,  0, 0,  2, -1, 3, -1,
5959371c9d4SSatish Balay                                               /*  1: b */};
596b5a892a1SMatthew G. Knepley   /* https://en.wikiversity.org/wiki/Symmetric_group_S4 */
5979371c9d4SSatish Balay   static const PetscInt tetArr[24 * 4 * 2] = {3, -2, 2, -3, 0, -1, 1, -1, /* -12: (1324)   p22 */
598b5a892a1SMatthew G. Knepley                                               3, -1, 1, -3, 2, -1, 0, -1, /* -11: (14)     p21 */
599b5a892a1SMatthew G. Knepley                                               3, -3, 0, -3, 1, -1, 2, -1, /* -10: (1234)   p18 */
600b5a892a1SMatthew G. Knepley                                               2, -1, 3, -1, 1, -3, 0, -2, /*  -9: (1423)   p17 */
601b5a892a1SMatthew G. Knepley                                               2, -3, 0, -1, 3, -2, 1, -3, /*  -8: (1342)   p13 */
602b5a892a1SMatthew G. Knepley                                               2, -2, 1, -2, 0, -2, 3, -2, /*  -7: (24)     p14 */
603b5a892a1SMatthew G. Knepley                                               1, -2, 0, -2, 2, -2, 3, -1, /*  -6: (34)     p6  */
604b5a892a1SMatthew G. Knepley                                               1, -1, 3, -3, 0, -3, 2, -2, /*  -5: (1243)   p10 */
605b5a892a1SMatthew G. Knepley                                               1, -3, 2, -1, 3, -1, 0, -3, /*  -4: (1432)   p9  */
606b5a892a1SMatthew G. Knepley                                               0, -3, 1, -1, 3, -3, 2, -3, /*  -3: (12)     p1  */
607b5a892a1SMatthew G. Knepley                                               0, -2, 2, -2, 1, -2, 3, -3, /*  -2: (23)     p2  */
608b5a892a1SMatthew G. Knepley                                               0, -1, 3, -2, 2, -3, 1, -2, /*  -1: (13)     p5  */
609b5a892a1SMatthew G. Knepley                                               0, 0,  1, 0,  2, 0,  3, 0,  /*   0: ()       p0  */
610b5a892a1SMatthew G. Knepley                                               0, 1,  3, 1,  1, 2,  2, 0,  /*   1: (123)    p4  */
611b5a892a1SMatthew G. Knepley                                               0, 2,  2, 1,  3, 0,  1, 2,  /*   2: (132)    p3  */
612b5a892a1SMatthew G. Knepley                                               1, 2,  0, 1,  3, 1,  2, 2,  /*   3: (12)(34) p7  */
613b5a892a1SMatthew G. Knepley                                               1, 0,  2, 0,  0, 0,  3, 1,  /*   4: (243)    p8  */
614b5a892a1SMatthew G. Knepley                                               1, 1,  3, 2,  2, 2,  0, 0,  /*   5: (143)    p11 */
615b5a892a1SMatthew G. Knepley                                               2, 1,  3, 0,  0, 2,  1, 0,  /*   6: (13)(24) p16 */
616b5a892a1SMatthew G. Knepley                                               2, 2,  1, 1,  3, 2,  0, 2,  /*   7: (142)    p15 */
617b5a892a1SMatthew G. Knepley                                               2, 0,  0, 0,  1, 0,  3, 2,  /*   8: (234)    p12 */
618b5a892a1SMatthew G. Knepley                                               3, 2,  2, 2,  1, 1,  0, 1,  /*   9: (14)(23) p23 */
619b5a892a1SMatthew G. Knepley                                               3, 0,  0, 2,  2, 1,  1, 1,  /*  10: (134)    p19 */
620b5a892a1SMatthew G. Knepley                                               3, 1,  1, 2,  0, 1,  2, 1 /*  11: (124)    p20 */};
621b5a892a1SMatthew G. Knepley   /* Each rotation determines a permutation of the four diagonals, and this defines the isomorphism with S_4 */
622ef8b56bfSJed Brown   static const PetscInt hexArr[48 * 6 * 2] = {
623b5a892a1SMatthew G. Knepley     2, -3, 3, -2, 4, -2, 5, -3, 1, -3, 0, -1, /* -24: reflect bottom and use -3 on top */
624b5a892a1SMatthew G. Knepley     4, -2, 5, -2, 0, -1, 1, -4, 3, -2, 2, -3, /* -23: reflect bottom and use -3 on top */
625b5a892a1SMatthew G. Knepley     5, -3, 4, -1, 1, -2, 0, -3, 3, -4, 2, -1, /* -22: reflect bottom and use -3 on top */
626b5a892a1SMatthew G. Knepley     3, -1, 2, -4, 4, -4, 5, -1, 0, -4, 1, -4, /* -21: reflect bottom and use -3 on top */
627b5a892a1SMatthew G. Knepley     3, -3, 2, -2, 5, -1, 4, -4, 1, -1, 0, -3, /* -20: reflect bottom and use -3 on top */
628b5a892a1SMatthew G. Knepley     4, -4, 5, -4, 1, -4, 0, -1, 2, -4, 3, -1, /* -19: reflect bottom and use -3 on top */
629b5a892a1SMatthew G. Knepley     2, -1, 3, -4, 5, -3, 4, -2, 0, -2, 1, -2, /* -18: reflect bottom and use -3 on top */
630b5a892a1SMatthew G. Knepley     5, -1, 4, -3, 0, -3, 1, -2, 2, -2, 3, -3, /* -17: reflect bottom and use -3 on top */
631b5a892a1SMatthew G. Knepley     4, -3, 5, -1, 3, -2, 2, -4, 1, -4, 0, -4, /* -16: reflect bottom and use -3 on top */
632b5a892a1SMatthew G. Knepley     5, -4, 4, -4, 3, -4, 2, -2, 0, -3, 1, -1, /* -15: reflect bottom and use -3 on top */
633b5a892a1SMatthew G. Knepley     3, -4, 2, -1, 1, -1, 0, -4, 4, -4, 5, -4, /* -14: reflect bottom and use -3 on top */
634b5a892a1SMatthew G. Knepley     2, -2, 3, -3, 0, -2, 1, -3, 4, -2, 5, -2, /* -13: reflect bottom and use -3 on top */
635b5a892a1SMatthew G. Knepley     1, -3, 0, -1, 4, -1, 5, -4, 3, -1, 2, -4, /* -12: reflect bottom and use -3 on top */
636b5a892a1SMatthew G. Knepley     1, -1, 0, -3, 5, -4, 4, -1, 2, -1, 3, -4, /* -11: reflect bottom and use -3 on top */
637b5a892a1SMatthew G. Knepley     5, -2, 4, -2, 2, -2, 3, -4, 1, -2, 0, -2, /* -10: reflect bottom and use -3 on top */
638b5a892a1SMatthew G. Knepley     1, -2, 0, -2, 2, -1, 3, -1, 4, -1, 5, -3, /*  -9: reflect bottom and use -3 on top */
639b5a892a1SMatthew G. Knepley     4, -1, 5, -3, 2, -4, 3, -2, 0, -1, 1, -3, /*  -8: reflect bottom and use -3 on top */
640b5a892a1SMatthew G. Knepley     3, -2, 2, -3, 0, -4, 1, -1, 5, -1, 4, -3, /*  -7: reflect bottom and use -3 on top */
641b5a892a1SMatthew G. Knepley     1, -4, 0, -4, 3, -1, 2, -1, 5, -4, 4, -4, /*  -6: reflect bottom and use -3 on top */
642b5a892a1SMatthew G. Knepley     2, -4, 3, -1, 1, -3, 0, -2, 5, -3, 4, -1, /*  -5: reflect bottom and use -3 on top */
643b5a892a1SMatthew G. Knepley     0, -4, 1, -4, 4, -3, 5, -2, 2, -3, 3, -2, /*  -4: reflect bottom and use -3 on top */
644b5a892a1SMatthew G. Knepley     0, -3, 1, -1, 3, -3, 2, -3, 4, -3, 5, -1, /*  -3: reflect bottom and use -3 on top */
645b5a892a1SMatthew G. Knepley     0, -2, 1, -2, 5, -2, 4, -3, 3, -3, 2, -2, /*  -2: reflect bottom and use -3 on top */
646b5a892a1SMatthew G. Knepley     0, -1, 1, -3, 2, -3, 3, -3, 5, -2, 4, -2, /*  -1: reflect bottom and use -3 on top */
647b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  5, 0,  /*   0: identity */
648b5a892a1SMatthew G. Knepley     0, 1,  1, 3,  5, 3,  4, 0,  2, 0,  3, 1,  /*   1: 90  rotation about z */
649b5a892a1SMatthew G. Knepley     0, 2,  1, 2,  3, 0,  2, 0,  5, 3,  4, 1,  /*   2: 180 rotation about z */
650b5a892a1SMatthew G. Knepley     0, 3,  1, 1,  4, 0,  5, 3,  3, 0,  2, 1,  /*   3: 270 rotation about z */
651b5a892a1SMatthew G. Knepley     2, 3,  3, 2,  1, 0,  0, 3,  4, 3,  5, 1,  /*   4: 90  rotation about x */
652b5a892a1SMatthew G. Knepley     1, 3,  0, 1,  3, 2,  2, 2,  4, 2,  5, 2,  /*   5: 180 rotation about x */
653b5a892a1SMatthew G. Knepley     3, 1,  2, 0,  0, 1,  1, 2,  4, 1,  5, 3,  /*   6: 270 rotation about x */
654b5a892a1SMatthew G. Knepley     4, 0,  5, 0,  2, 1,  3, 3,  1, 1,  0, 3,  /*   7: 90  rotation about y */
655b5a892a1SMatthew G. Knepley     1, 1,  0, 3,  2, 2,  3, 2,  5, 1,  4, 3,  /*   8: 180 rotation about y */
656b5a892a1SMatthew G. Knepley     5, 1,  4, 3,  2, 3,  3, 1,  0, 0,  1, 0,  /*   9: 270 rotation about y */
657b5a892a1SMatthew G. Knepley     1, 0,  0, 0,  5, 1,  4, 2,  3, 2,  2, 3,  /*  10: 180 rotation about x+y */
658b5a892a1SMatthew G. Knepley     1, 2,  0, 2,  4, 2,  5, 1,  2, 2,  3, 3,  /*  11: 180 rotation about x-y */
659b5a892a1SMatthew G. Knepley     2, 1,  3, 0,  0, 3,  1, 0,  5, 0,  4, 0,  /*  12: 180 rotation about y+z */
660b5a892a1SMatthew G. Knepley     3, 3,  2, 2,  1, 2,  0, 1,  5, 2,  4, 2,  /*  13: 180 rotation about y-z */
661b5a892a1SMatthew G. Knepley     5, 3,  4, 1,  3, 1,  2, 3,  1, 3,  0, 1,  /*  14: 180 rotation about z+x */
662b5a892a1SMatthew G. Knepley     4, 2,  5, 2,  3, 3,  2, 1,  0, 2,  1, 2,  /*  15: 180 rotation about z-x */
663b5a892a1SMatthew G. Knepley     5, 0,  4, 0,  0, 0,  1, 3,  3, 1,  2, 0,  /*  16: 120 rotation about x+y+z (v0v6) */
664b5a892a1SMatthew G. Knepley     2, 0,  3, 1,  5, 0,  4, 3,  1, 0,  0, 0,  /*  17: 240 rotation about x+y+z (v0v6) */
665b5a892a1SMatthew G. Knepley     4, 3,  5, 1,  1, 1,  0, 2,  3, 3,  2, 2,  /*  18: 120 rotation about x+y-z (v4v2) */
666b5a892a1SMatthew G. Knepley     3, 2,  2, 3,  5, 2,  4, 1,  0, 1,  1, 3,  /*  19: 240 rotation about x+y-z (v4v2) */
667b5a892a1SMatthew G. Knepley     3, 0,  2, 1,  4, 1,  5, 2,  1, 2,  0, 2,  /*  20: 120 rotation about x-y+z (v1v5) */
668b5a892a1SMatthew G. Knepley     5, 2,  4, 2,  1, 3,  0, 0,  2, 3,  3, 2,  /*  21: 240 rotation about x-y+z (v1v5) */
669b5a892a1SMatthew G. Knepley     4, 1,  5, 3,  0, 2,  1, 1,  2, 1,  3, 0,  /*  22: 120 rotation about x-y-z (v7v3) */
670b5a892a1SMatthew G. Knepley     2, 2,  3, 3,  4, 3,  5, 0,  0, 3,  1, 1,  /*  23: 240 rotation about x-y-z (v7v3) */
671b5a892a1SMatthew G. Knepley   };
672ef8b56bfSJed Brown   static const PetscInt tripArr[12 * 5 * 2] = {
673b5a892a1SMatthew G. Knepley     1, -3, 0, -1, 3, -1, 4, -1, 2, -1, /* -6: reflect bottom and top */
674b5a892a1SMatthew G. Knepley     1, -1, 0, -3, 4, -1, 2, -1, 3, -1, /* -5: reflect bottom and top */
675b5a892a1SMatthew G. Knepley     1, -2, 0, -2, 2, -1, 3, -1, 4, -1, /* -4: reflect bottom and top */
676b5a892a1SMatthew G. Knepley     0, -3, 1, -1, 3, -3, 2, -3, 4, -3, /* -3: reflect bottom and top */
677b5a892a1SMatthew G. Knepley     0, -2, 1, -2, 4, -3, 3, -3, 2, -3, /* -2: reflect bottom and top */
678b5a892a1SMatthew G. Knepley     0, -1, 1, -3, 2, -3, 4, -3, 3, -3, /* -1: reflect bottom and top */
679b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  /*  0: identity */
680b5a892a1SMatthew G. Knepley     0, 1,  1, 2,  4, 0,  2, 0,  3, 0,  /*  1: 120 rotation about z */
681b5a892a1SMatthew G. Knepley     0, 2,  1, 1,  3, 0,  4, 0,  2, 0,  /*  2: 240 rotation about z */
682b5a892a1SMatthew G. Knepley     1, 1,  0, 2,  2, 2,  4, 2,  3, 2,  /*  3: 180 rotation about y of 0 */
683b5a892a1SMatthew G. Knepley     1, 0,  0, 0,  4, 2,  3, 2,  2, 2,  /*  4: 180 rotation about y of 1 */
684b5a892a1SMatthew G. Knepley     1, 2,  0, 1,  3, 2,  2, 2,  4, 2,  /*  5: 180 rotation about y of 2 */
685b5a892a1SMatthew G. Knepley   };
686b5a892a1SMatthew G. Knepley   /* a: rotate 120 about z
687b5a892a1SMatthew G. Knepley      b: swap top and bottom segments
688b5a892a1SMatthew G. Knepley      r: reflect */
689ef8b56bfSJed Brown   static const PetscInt ttriArr[12 * 5 * 2] = {
690b5a892a1SMatthew G. Knepley     1, -3, 0, -3, 2, -2, 4, -2, 3, -2, /* -6: r b a^2 */
691b5a892a1SMatthew G. Knepley     1, -2, 0, -2, 4, -2, 3, -2, 2, -2, /* -5: r b a */
692b5a892a1SMatthew G. Knepley     1, -1, 0, -1, 3, -2, 2, -2, 4, -2, /* -4: r b */
693b5a892a1SMatthew G. Knepley     0, -3, 1, -3, 2, -1, 4, -1, 3, -1, /* -3: r a^2 */
694b5a892a1SMatthew G. Knepley     0, -2, 1, -2, 4, -1, 3, -1, 2, -1, /* -2: r a */
695b5a892a1SMatthew G. Knepley     0, -1, 1, -1, 3, -1, 2, -1, 4, -1, /* -1: r */
696b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  /*  0: identity */
697b5a892a1SMatthew G. Knepley     0, 1,  1, 1,  3, 0,  4, 0,  2, 0,  /*  1: a */
698b5a892a1SMatthew G. Knepley     0, 2,  1, 2,  4, 0,  2, 0,  3, 0,  /*  2: a^2 */
699b5a892a1SMatthew G. Knepley     1, 0,  0, 0,  2, 1,  3, 1,  4, 1,  /*  3: b */
700b5a892a1SMatthew G. Knepley     1, 1,  0, 1,  3, 1,  4, 1,  2, 1,  /*  4: b a */
701b5a892a1SMatthew G. Knepley     1, 2,  0, 2,  4, 1,  2, 1,  3, 1,  /*  5: b a^2 */
702b5a892a1SMatthew G. Knepley   };
703b5a892a1SMatthew G. Knepley   /* a: rotate 90 about z
704b5a892a1SMatthew G. Knepley      b: swap top and bottom segments
705b5a892a1SMatthew G. Knepley      r: reflect */
706ef8b56bfSJed Brown   static const PetscInt tquadArr[16 * 6 * 2] = {
707b5a892a1SMatthew G. Knepley     1, -4, 0, -4, 3, -2, 2, -2, 5, -2, 4, -2, /* -8: r b a^3 */
708b5a892a1SMatthew G. Knepley     1, -3, 0, -3, 2, -2, 5, -2, 4, -2, 3, -2, /* -7: r b a^2 */
709b5a892a1SMatthew G. Knepley     1, -2, 0, -2, 5, -2, 4, -2, 3, -2, 2, -2, /* -6: r b a */
710b5a892a1SMatthew G. Knepley     1, -1, 0, -1, 4, -2, 3, -2, 2, -2, 5, -2, /* -5: r b */
711b5a892a1SMatthew G. Knepley     0, -4, 1, -4, 3, -1, 2, -1, 5, -1, 4, -1, /* -4: r a^3 */
712b5a892a1SMatthew G. Knepley     0, -3, 1, -3, 2, -1, 5, -1, 4, -1, 3, -1, /* -3: r a^2 */
713b5a892a1SMatthew G. Knepley     0, -2, 1, -2, 5, -1, 4, -1, 3, -1, 2, -1, /* -2: r a */
714b5a892a1SMatthew G. Knepley     0, -1, 1, -1, 4, -1, 3, -1, 2, -1, 5, -1, /* -1: r */
715b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  5, 0,  /*  0: identity */
716b5a892a1SMatthew G. Knepley     0, 1,  1, 1,  3, 0,  4, 0,  5, 0,  2, 0,  /*  1: a */
717b5a892a1SMatthew G. Knepley     0, 2,  1, 2,  4, 0,  5, 0,  2, 0,  3, 0,  /*  2: a^2 */
718b5a892a1SMatthew G. Knepley     0, 3,  1, 3,  5, 0,  2, 0,  3, 0,  4, 0,  /*  3: a^3 */
719b5a892a1SMatthew G. Knepley     1, 0,  0, 0,  2, 1,  3, 1,  4, 1,  5, 1,  /*  4: b */
720b5a892a1SMatthew G. Knepley     1, 1,  0, 1,  3, 1,  4, 1,  5, 1,  2, 1,  /*  5: b a */
721b5a892a1SMatthew G. Knepley     1, 2,  0, 2,  4, 1,  5, 1,  2, 1,  3, 1,  /*  6: b a^2 */
722b5a892a1SMatthew G. Knepley     1, 3,  0, 3,  5, 1,  2, 1,  3, 1,  4, 1,  /*  7: b a^3 */
723b5a892a1SMatthew G. Knepley   };
724ef8b56bfSJed Brown   static const PetscInt pyrArr[8 * 5 * 2] = {
725b5a892a1SMatthew G. Knepley     0, -4, 2, -3, 1, -3, 4, -3, 3, -3, /* -4: Reflect bottom face */
726b5a892a1SMatthew G. Knepley     0, -3, 3, -3, 2, -3, 1, -3, 4, -3, /* -3: Reflect bottom face */
727b5a892a1SMatthew G. Knepley     0, -2, 4, -3, 3, -3, 2, -3, 1, -3, /* -2: Reflect bottom face */
728b5a892a1SMatthew G. Knepley     0, -1, 1, -3, 4, -3, 3, -3, 2, -3, /* -1: Reflect bottom face */
729b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  /*  0: identity */
730b5a892a1SMatthew G. Knepley     0, 1,  4, 0,  1, 0,  2, 0,  3, 0,  /*  1:  90 rotation about z */
731b5a892a1SMatthew G. Knepley     0, 2,  3, 0,  4, 0,  1, 0,  2, 0,  /*  2: 180 rotation about z */
732b5a892a1SMatthew G. Knepley     0, 3,  2, 0,  3, 0,  4, 0,  1, 0,  /*  3: 270 rotation about z */
733b5a892a1SMatthew G. Knepley   };
734b5a892a1SMatthew G. Knepley   switch (ct) {
735d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
736d71ae5a4SJacob Faibussowitsch     return pntArr;
737d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
738d71ae5a4SJacob Faibussowitsch     return &segArr[(o + 1) * 2 * 2];
739d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
740d71ae5a4SJacob Faibussowitsch     return &segArr[(o + 1) * 2 * 2];
741d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
742d71ae5a4SJacob Faibussowitsch     return &triArr[(o + 3) * 3 * 2];
743d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
744d71ae5a4SJacob Faibussowitsch     return &quadArr[(o + 4) * 4 * 2];
745d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
746d71ae5a4SJacob Faibussowitsch     return &tsegArr[(o + 2) * 4 * 2];
747d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
748d71ae5a4SJacob Faibussowitsch     return &tetArr[(o + 12) * 4 * 2];
749d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
750d71ae5a4SJacob Faibussowitsch     return &hexArr[(o + 24) * 6 * 2];
751d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
752d71ae5a4SJacob Faibussowitsch     return &tripArr[(o + 6) * 5 * 2];
753d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
754d71ae5a4SJacob Faibussowitsch     return &ttriArr[(o + 6) * 5 * 2];
755d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
756d71ae5a4SJacob Faibussowitsch     return &tquadArr[(o + 8) * 6 * 2];
757d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
758d71ae5a4SJacob Faibussowitsch     return &pyrArr[(o + 4) * 5 * 2];
759d71ae5a4SJacob Faibussowitsch   default:
760f22e26b7SPierre Jolivet     return PETSC_NULLPTR;
761b5a892a1SMatthew G. Knepley   }
762b5a892a1SMatthew G. Knepley }
763b5a892a1SMatthew G. Knepley 
764d5b43468SJose E. Roman /* A vertex arrangement is a vertex order */
765d71ae5a4SJacob Faibussowitsch static inline const PetscInt *DMPolytopeTypeGetVertexArrangment(DMPolytopeType ct, PetscInt o)
766d71ae5a4SJacob Faibussowitsch {
767ef8b56bfSJed Brown   static const PetscInt pntVerts[1]      = {0};
7689371c9d4SSatish Balay   static const PetscInt segVerts[2 * 2]  = {1, 0, 0, 1};
7699371c9d4SSatish 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};
7709371c9d4SSatish 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};
7719371c9d4SSatish Balay   static const PetscInt tsegVerts[4 * 4] = {3, 2, 1, 0, 1, 0, 3, 2, 0, 1, 2, 3, 2, 3, 0, 1};
7729371c9d4SSatish Balay   static const PetscInt tetVerts[24 * 4] = {2, 3, 1, 0, /* -12: (1324)   p22 */
773b5a892a1SMatthew G. Knepley                                             3, 1, 2, 0, /* -11: (14)     p21 */
774b5a892a1SMatthew G. Knepley                                             1, 2, 3, 0, /* -10: (1234)   p18 */
775b5a892a1SMatthew G. Knepley                                             3, 2, 0, 1, /*  -9: (1423)   p17 */
776b5a892a1SMatthew G. Knepley                                             2, 0, 3, 1, /*  -8: (1342)   p13 */
777b5a892a1SMatthew G. Knepley                                             0, 3, 2, 1, /*  -7: (24)     p14 */
778b5a892a1SMatthew G. Knepley                                             0, 1, 3, 2, /*  -6: (34)     p6  */
779b5a892a1SMatthew G. Knepley                                             1, 3, 0, 2, /*  -5: (1243)   p10 */
780b5a892a1SMatthew G. Knepley                                             3, 0, 1, 2, /*  -4: (1432    p9  */
781b5a892a1SMatthew G. Knepley                                             1, 0, 2, 3, /*  -3: (12)     p1  */
782b5a892a1SMatthew G. Knepley                                             0, 2, 1, 3, /*  -2: (23)     p2  */
783b5a892a1SMatthew G. Knepley                                             2, 1, 0, 3, /*  -1: (13)     p5  */
784b5a892a1SMatthew G. Knepley                                             0, 1, 2, 3, /*   0: ()       p0  */
785b5a892a1SMatthew G. Knepley                                             1, 2, 0, 3, /*   1: (123)    p4  */
786b5a892a1SMatthew G. Knepley                                             2, 0, 1, 3, /*   2: (132)    p3  */
787b5a892a1SMatthew G. Knepley                                             1, 0, 3, 2, /*   3: (12)(34) p7  */
788b5a892a1SMatthew G. Knepley                                             0, 3, 1, 2, /*   4: (243)    p8  */
789b5a892a1SMatthew G. Knepley                                             3, 1, 0, 2, /*   5: (143)    p11 */
790b5a892a1SMatthew G. Knepley                                             2, 3, 0, 1, /*   6: (13)(24) p16 */
791b5a892a1SMatthew G. Knepley                                             3, 0, 2, 1, /*   7: (142)    p15 */
792b5a892a1SMatthew G. Knepley                                             0, 2, 3, 1, /*   8: (234)    p12 */
793b5a892a1SMatthew G. Knepley                                             3, 2, 1, 0, /*   9: (14)(23) p23 */
794b5a892a1SMatthew G. Knepley                                             2, 1, 3, 0, /*  10: (134)    p19 */
795b5a892a1SMatthew G. Knepley                                             1, 3, 2, 0 /*  11: (124)    p20 */};
796ef8b56bfSJed Brown   static const PetscInt hexVerts[48 * 8] = {
797b5a892a1SMatthew G. Knepley     3, 0, 4, 5, 2, 6, 7, 1, /* -24: reflected 23 */
798b5a892a1SMatthew G. Knepley     3, 5, 6, 2, 0, 1, 7, 4, /* -23: reflected 22 */
799b5a892a1SMatthew G. Knepley     4, 0, 1, 7, 5, 6, 2, 3, /* -22: reflected 21 */
800b5a892a1SMatthew G. Knepley     6, 7, 1, 2, 5, 3, 0, 4, /* -21: reflected 20 */
801b5a892a1SMatthew G. Knepley     1, 2, 6, 7, 0, 4, 5, 3, /* -20: reflected 19 */
802b5a892a1SMatthew G. Knepley     6, 2, 3, 5, 7, 4, 0, 1, /* -19: reflected 18 */
803b5a892a1SMatthew G. Knepley     4, 5, 3, 0, 7, 1, 2, 6, /* -18: reflected 17 */
804b5a892a1SMatthew G. Knepley     1, 7, 4, 0, 2, 3, 5, 6, /* -17: reflected 16 */
805b5a892a1SMatthew G. Knepley     2, 3, 5, 6, 1, 7, 4, 0, /* -16: reflected 15 */
806b5a892a1SMatthew G. Knepley     7, 4, 0, 1, 6, 2, 3, 5, /* -15: reflected 14 */
807b5a892a1SMatthew G. Knepley     7, 1, 2, 6, 4, 5, 3, 0, /* -14: reflected 13 */
808b5a892a1SMatthew G. Knepley     0, 4, 5, 3, 1, 2, 6, 7, /* -13: reflected 12 */
809b5a892a1SMatthew G. Knepley     5, 4, 7, 6, 3, 2, 1, 0, /* -12: reflected 11 */
810b5a892a1SMatthew G. Knepley     7, 6, 5, 4, 1, 0, 3, 2, /* -11: reflected 10 */
811b5a892a1SMatthew G. Knepley     0, 1, 7, 4, 3, 5, 6, 2, /* -10: reflected  9 */
812b5a892a1SMatthew G. Knepley     4, 7, 6, 5, 0, 3, 2, 1, /*  -9: reflected  8 */
813b5a892a1SMatthew G. Knepley     5, 6, 2, 3, 4, 0, 1, 7, /*  -8: reflected  7 */
814b5a892a1SMatthew G. Knepley     2, 6, 7, 1, 3, 0, 4, 5, /*  -7: reflected  6 */
815b5a892a1SMatthew G. Knepley     6, 5, 4, 7, 2, 1, 0, 3, /*  -6: reflected  5 */
816b5a892a1SMatthew G. Knepley     5, 3, 0, 4, 6, 7, 1, 2, /*  -5: reflected  4 */
817b5a892a1SMatthew G. Knepley     2, 1, 0, 3, 6, 5, 4, 7, /*  -4: reflected  3 */
818b5a892a1SMatthew G. Knepley     1, 0, 3, 2, 7, 6, 5, 4, /*  -3: reflected  2 */
819b5a892a1SMatthew G. Knepley     0, 3, 2, 1, 4, 7, 6, 5, /*  -2: reflected  1 */
820b5a892a1SMatthew G. Knepley     3, 2, 1, 0, 5, 4, 7, 6, /*  -1: reflected  0 */
821b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, 5, 6, 7, /*   0: identity */
822b5a892a1SMatthew G. Knepley     1, 2, 3, 0, 7, 4, 5, 6, /*   1: 90  rotation about z */
823b5a892a1SMatthew G. Knepley     2, 3, 0, 1, 6, 7, 4, 5, /*   2: 180 rotation about z */
824b5a892a1SMatthew G. Knepley     3, 0, 1, 2, 5, 6, 7, 4, /*   3: 270 rotation about z */
825b5a892a1SMatthew G. Knepley     4, 0, 3, 5, 7, 6, 2, 1, /*   4: 90  rotation about x */
826b5a892a1SMatthew G. Knepley     7, 4, 5, 6, 1, 2, 3, 0, /*   5: 180 rotation about x */
827b5a892a1SMatthew G. Knepley     1, 7, 6, 2, 0, 3, 5, 4, /*   6: 270 rotation about x */
828b5a892a1SMatthew G. Knepley     3, 2, 6, 5, 0, 4, 7, 1, /*   7: 90  rotation about y */
829b5a892a1SMatthew G. Knepley     5, 6, 7, 4, 3, 0, 1, 2, /*   8: 180 rotation about y */
830b5a892a1SMatthew G. Knepley     4, 7, 1, 0, 5, 3, 2, 6, /*   9: 270 rotation about y */
831b5a892a1SMatthew G. Knepley     4, 5, 6, 7, 0, 1, 2, 3, /*  10: 180 rotation about x+y */
832b5a892a1SMatthew G. Knepley     6, 7, 4, 5, 2, 3, 0, 1, /*  11: 180 rotation about x-y */
833b5a892a1SMatthew G. Knepley     3, 5, 4, 0, 2, 1, 7, 6, /*  12: 180 rotation about y+z */
834b5a892a1SMatthew G. Knepley     6, 2, 1, 7, 5, 4, 0, 3, /*  13: 180 rotation about y-z */
835b5a892a1SMatthew G. Knepley     1, 0, 4, 7, 2, 6, 5, 3, /*  14: 180 rotation about z+x */
836b5a892a1SMatthew G. Knepley     6, 5, 3, 2, 7, 1, 0, 4, /*  15: 180 rotation about z-x */
837b5a892a1SMatthew G. Knepley     0, 4, 7, 1, 3, 2, 6, 5, /*  16: 120 rotation about x+y+z (v0v6) */
838b5a892a1SMatthew G. Knepley     0, 3, 5, 4, 1, 7, 6, 2, /*  17: 240 rotation about x+y+z (v0v6) */
839b5a892a1SMatthew G. Knepley     5, 3, 2, 6, 4, 7, 1, 0, /*  18: 120 rotation about x+y-z (v4v2) */
840b5a892a1SMatthew G. Knepley     7, 6, 2, 1, 4, 0, 3, 5, /*  19: 240 rotation about x+y-z (v4v2) */
841b5a892a1SMatthew G. Knepley     2, 1, 7, 6, 3, 5, 4, 0, /*  20: 120 rotation about x-y+z (v1v5) */
842b5a892a1SMatthew G. Knepley     7, 1, 0, 4, 6, 5, 3, 2, /*  21: 240 rotation about x-y+z (v1v5) */
843b5a892a1SMatthew G. Knepley     2, 6, 5, 3, 1, 0, 4, 7, /*  22: 120 rotation about x-y-z (v7v3) */
844b5a892a1SMatthew G. Knepley     5, 4, 0, 3, 6, 2, 1, 7, /*  23: 240 rotation about x-y-z (v7v3) */
845b5a892a1SMatthew G. Knepley   };
846ef8b56bfSJed Brown   static const PetscInt tripVerts[12 * 6] = {
847b5a892a1SMatthew G. Knepley     4, 3, 5, 2, 1, 0, /* -6: reflect bottom and top */
848b5a892a1SMatthew G. Knepley     5, 4, 3, 1, 0, 2, /* -5: reflect bottom and top */
849b5a892a1SMatthew G. Knepley     3, 5, 4, 0, 2, 1, /* -4: reflect bottom and top */
850b5a892a1SMatthew G. Knepley     1, 0, 2, 5, 4, 3, /* -3: reflect bottom and top */
851b5a892a1SMatthew G. Knepley     0, 2, 1, 3, 5, 4, /* -2: reflect bottom and top */
852b5a892a1SMatthew G. Knepley     2, 1, 0, 4, 3, 5, /* -1: reflect bottom and top */
853b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, 5, /*  0: identity */
854b5a892a1SMatthew G. Knepley     1, 2, 0, 5, 3, 4, /*  1: 120 rotation about z */
855b5a892a1SMatthew G. Knepley     2, 0, 1, 4, 5, 3, /*  2: 240 rotation about z */
856b5a892a1SMatthew G. Knepley     4, 5, 3, 2, 0, 1, /*  3: 180 rotation about y of 0 */
857b5a892a1SMatthew G. Knepley     3, 4, 5, 0, 1, 2, /*  4: 180 rotation about y of 1 */
858b5a892a1SMatthew G. Knepley     5, 3, 4, 1, 2, 0, /*  5: 180 rotation about y of 2 */
859b5a892a1SMatthew G. Knepley   };
860ef8b56bfSJed Brown   static const PetscInt ttriVerts[12 * 6] = {
861b5a892a1SMatthew G. Knepley     4, 3, 5, 1, 0, 2, /* -6: r b a^2 */
862b5a892a1SMatthew G. Knepley     3, 5, 4, 0, 2, 1, /* -5: r b a */
863b5a892a1SMatthew G. Knepley     5, 4, 3, 2, 1, 0, /* -4: r b */
864b5a892a1SMatthew G. Knepley     1, 0, 2, 4, 3, 5, /* -3: r a^2 */
865b5a892a1SMatthew G. Knepley     0, 2, 1, 3, 5, 4, /* -2: r a */
866b5a892a1SMatthew G. Knepley     2, 1, 0, 5, 4, 3, /* -1: r */
867b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, 5, /*  0: identity */
868b5a892a1SMatthew G. Knepley     1, 2, 0, 4, 5, 3, /*  1: a */
869b5a892a1SMatthew G. Knepley     2, 0, 1, 5, 3, 4, /*  2: a^2 */
870b5a892a1SMatthew G. Knepley     3, 4, 5, 0, 1, 2, /*  3: b */
871b5a892a1SMatthew G. Knepley     4, 5, 3, 1, 2, 0, /*  4: b a */
872b5a892a1SMatthew G. Knepley     5, 3, 4, 2, 0, 1, /*  5: b a^2 */
873b5a892a1SMatthew G. Knepley   };
874b5a892a1SMatthew G. Knepley   /* a: rotate 90 about z
875b5a892a1SMatthew G. Knepley      b: swap top and bottom segments
876b5a892a1SMatthew G. Knepley      r: reflect */
877ef8b56bfSJed Brown   static const PetscInt tquadVerts[16 * 8] = {
878b5a892a1SMatthew G. Knepley     6, 5, 4, 7, 2, 1, 0, 3, /* -8: r b a^3 */
879b5a892a1SMatthew G. Knepley     5, 4, 7, 6, 1, 0, 3, 2, /* -7: r b a^2 */
880b5a892a1SMatthew G. Knepley     4, 7, 6, 5, 0, 3, 2, 1, /* -6: r b a */
881b5a892a1SMatthew G. Knepley     7, 6, 5, 4, 3, 2, 1, 0, /* -5: r b */
882b5a892a1SMatthew G. Knepley     2, 1, 0, 3, 6, 5, 4, 7, /* -4: r a^3 */
883b5a892a1SMatthew G. Knepley     1, 0, 3, 2, 5, 4, 7, 6, /* -3: r a^2 */
884b5a892a1SMatthew G. Knepley     0, 3, 2, 1, 4, 7, 6, 5, /* -2: r a */
885b5a892a1SMatthew G. Knepley     3, 2, 1, 0, 7, 6, 5, 4, /* -1: r */
886b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, 5, 6, 7, /*  0: identity */
887b5a892a1SMatthew G. Knepley     1, 2, 3, 0, 5, 6, 7, 4, /*  1: a */
888b5a892a1SMatthew G. Knepley     2, 3, 0, 1, 6, 7, 4, 5, /*  2: a^2 */
889b5a892a1SMatthew G. Knepley     3, 0, 1, 2, 7, 4, 5, 6, /*  3: a^3 */
890b5a892a1SMatthew G. Knepley     4, 5, 6, 7, 0, 1, 2, 3, /*  4: b */
891b5a892a1SMatthew G. Knepley     5, 6, 7, 4, 1, 2, 3, 0, /*  5: b a */
892b5a892a1SMatthew G. Knepley     6, 7, 4, 5, 2, 3, 0, 1, /*  6: b a^2 */
893b5a892a1SMatthew G. Knepley     7, 4, 5, 6, 3, 0, 1, 2, /*  7: b a^3 */
894b5a892a1SMatthew G. Knepley   };
895ef8b56bfSJed Brown   static const PetscInt pyrVerts[8 * 5] = {
896b5a892a1SMatthew G. Knepley     2, 1, 0, 3, 4, /* -4: Reflect bottom face */
897b5a892a1SMatthew G. Knepley     1, 0, 3, 2, 4, /* -3: Reflect bottom face */
898b5a892a1SMatthew G. Knepley     0, 3, 2, 1, 4, /* -2: Reflect bottom face */
899b5a892a1SMatthew G. Knepley     3, 2, 1, 0, 4, /* -1: Reflect bottom face */
900b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, /*  0: identity */
901b5a892a1SMatthew G. Knepley     1, 2, 3, 0, 4, /*  1:  90 rotation about z */
902b5a892a1SMatthew G. Knepley     2, 3, 0, 1, 4, /*  2: 180 rotation about z */
903b5a892a1SMatthew G. Knepley     3, 0, 1, 2, 4, /*  3: 270 rotation about z */
904b5a892a1SMatthew G. Knepley   };
905b5a892a1SMatthew G. Knepley   switch (ct) {
906d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
907d71ae5a4SJacob Faibussowitsch     return pntVerts;
908d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
909d71ae5a4SJacob Faibussowitsch     return &segVerts[(o + 1) * 2];
910d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
911d71ae5a4SJacob Faibussowitsch     return &segVerts[(o + 1) * 2];
912d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
913d71ae5a4SJacob Faibussowitsch     return &triVerts[(o + 3) * 3];
914d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
915d71ae5a4SJacob Faibussowitsch     return &quadVerts[(o + 4) * 4];
916d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
917d71ae5a4SJacob Faibussowitsch     return &tsegVerts[(o + 2) * 4];
918d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
919d71ae5a4SJacob Faibussowitsch     return &tetVerts[(o + 12) * 4];
920d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
921d71ae5a4SJacob Faibussowitsch     return &hexVerts[(o + 24) * 8];
922d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
923d71ae5a4SJacob Faibussowitsch     return &tripVerts[(o + 6) * 6];
924d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
925d71ae5a4SJacob Faibussowitsch     return &ttriVerts[(o + 6) * 6];
926d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
927d71ae5a4SJacob Faibussowitsch     return &tquadVerts[(o + 8) * 8];
928d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
929d71ae5a4SJacob Faibussowitsch     return &pyrVerts[(o + 4) * 5];
930d71ae5a4SJacob Faibussowitsch   default:
931f22e26b7SPierre Jolivet     return PETSC_NULLPTR;
932b5a892a1SMatthew G. Knepley   }
933b5a892a1SMatthew G. Knepley }
934b5a892a1SMatthew G. Knepley 
935b5a892a1SMatthew G. Knepley /* This is orientation o1 acting on orientation o2 */
936d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeComposeOrientation(DMPolytopeType ct, PetscInt o1, PetscInt o2)
937d71ae5a4SJacob Faibussowitsch {
9389371c9d4SSatish Balay   static const PetscInt segMult[2 * 2]   = {0, -1, -1, 0};
9399371c9d4SSatish 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};
9409371c9d4SSatish 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,
9419371c9d4SSatish 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};
9429371c9d4SSatish Balay   static const PetscInt tsegMult[4 * 4]  = {0, 1, -2, -1, 1, 0, -1, -2, -2, -1, 0, 1, -1, -2, 1, 0};
943ef8b56bfSJed Brown   static const PetscInt tetMult[24 * 24] = {
9449371c9d4SSatish 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,
9459371c9d4SSatish 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,
9469371c9d4SSatish 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,
9479371c9d4SSatish 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,
9489371c9d4SSatish 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,
9499371c9d4SSatish 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,
9509371c9d4SSatish 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,
9519371c9d4SSatish 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,
9529371c9d4SSatish 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,
9539371c9d4SSatish 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,
9549371c9d4SSatish 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,
9559371c9d4SSatish 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,
956b5a892a1SMatthew G. Knepley   };
957ef8b56bfSJed Brown   static const PetscInt hexMult[48 * 48] = {
958b5a892a1SMatthew 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,
959b5a892a1SMatthew 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,
960b5a892a1SMatthew 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,
961b5a892a1SMatthew 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,
962b5a892a1SMatthew 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,
963b5a892a1SMatthew 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,
964b5a892a1SMatthew 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,
965b5a892a1SMatthew 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,
966b5a892a1SMatthew 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,
967b5a892a1SMatthew 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,
968b5a892a1SMatthew 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,
969b5a892a1SMatthew 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,
970b5a892a1SMatthew 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,
971b5a892a1SMatthew 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,
972b5a892a1SMatthew 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,
973b5a892a1SMatthew 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,
974b5a892a1SMatthew 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,
975b5a892a1SMatthew 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,
976b5a892a1SMatthew 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,
977b5a892a1SMatthew 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,
978b5a892a1SMatthew 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,
979b5a892a1SMatthew 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,
980b5a892a1SMatthew 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,
981b5a892a1SMatthew 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,
982b5a892a1SMatthew 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,
983b5a892a1SMatthew 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,
984b5a892a1SMatthew 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,
985b5a892a1SMatthew 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,
986b5a892a1SMatthew 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,
987b5a892a1SMatthew 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,
988b5a892a1SMatthew 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,
989b5a892a1SMatthew 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,
990b5a892a1SMatthew 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,
991b5a892a1SMatthew 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,
992b5a892a1SMatthew 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,
993b5a892a1SMatthew 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,
994b5a892a1SMatthew 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,
995b5a892a1SMatthew 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,
996b5a892a1SMatthew 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,
997b5a892a1SMatthew 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,
998b5a892a1SMatthew 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,
999b5a892a1SMatthew 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,
1000b5a892a1SMatthew 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,
1001b5a892a1SMatthew 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,
1002b5a892a1SMatthew 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,
1003b5a892a1SMatthew 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,
1004b5a892a1SMatthew 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,
1005b5a892a1SMatthew 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,
1006b5a892a1SMatthew G. Knepley   };
1007ef8b56bfSJed Brown   static const PetscInt tripMult[12 * 12] = {
10089371c9d4SSatish 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,
10099371c9d4SSatish 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,
10109371c9d4SSatish 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,
1011b5a892a1SMatthew G. Knepley   };
1012ef8b56bfSJed Brown   static const PetscInt ttriMult[12 * 12] = {
10139371c9d4SSatish 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,
10149371c9d4SSatish 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,
10159371c9d4SSatish 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,
1016b5a892a1SMatthew G. Knepley   };
1017ef8b56bfSJed Brown   static const PetscInt tquadMult[16 * 16] = {
10189371c9d4SSatish 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,
10199371c9d4SSatish 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,
10209371c9d4SSatish 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,
10219371c9d4SSatish 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,
10229371c9d4SSatish 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,
1023b5a892a1SMatthew G. Knepley   };
1024ef8b56bfSJed Brown   static const PetscInt pyrMult[8 * 8] = {
10259371c9d4SSatish 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,
1026b5a892a1SMatthew G. Knepley   };
1027b5a892a1SMatthew G. Knepley   switch (ct) {
1028d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
1029d71ae5a4SJacob Faibussowitsch     return 0;
1030b5a892a1SMatthew G. Knepley   case DM_POLYTOPE_SEGMENT:
1031d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
1032d71ae5a4SJacob Faibussowitsch     return segMult[(o1 + 1) * 2 + o2 + 1];
1033d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
1034d71ae5a4SJacob Faibussowitsch     return triMult[(o1 + 3) * 6 + o2 + 3];
1035d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
1036d71ae5a4SJacob Faibussowitsch     return quadMult[(o1 + 4) * 8 + o2 + 4];
1037d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
1038d71ae5a4SJacob Faibussowitsch     return tsegMult[(o1 + 2) * 4 + o2 + 2];
1039d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
1040d71ae5a4SJacob Faibussowitsch     return tetMult[(o1 + 12) * 24 + o2 + 12];
1041d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
1042d71ae5a4SJacob Faibussowitsch     return hexMult[(o1 + 24) * 48 + o2 + 24];
1043d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
1044d71ae5a4SJacob Faibussowitsch     return tripMult[(o1 + 6) * 12 + o2 + 6];
1045d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
1046d71ae5a4SJacob Faibussowitsch     return ttriMult[(o1 + 6) * 12 + o2 + 6];
1047d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
1048d71ae5a4SJacob Faibussowitsch     return tquadMult[(o1 + 8) * 16 + o2 + 8];
1049d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
1050d71ae5a4SJacob Faibussowitsch     return pyrMult[(o1 + 4) * 8 + o2 + 4];
1051d71ae5a4SJacob Faibussowitsch   default:
1052d71ae5a4SJacob Faibussowitsch     return 0;
1053b5a892a1SMatthew G. Knepley   }
1054b5a892a1SMatthew G. Knepley }
1055b5a892a1SMatthew G. Knepley 
1056b5a892a1SMatthew G. Knepley /* This is orientation o1 acting on orientation o2^{-1} */
1057d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeComposeOrientationInv(DMPolytopeType ct, PetscInt o1, PetscInt o2)
1058d71ae5a4SJacob Faibussowitsch {
1059ef8b56bfSJed Brown   static const PetscInt triInv[6]    = {-3, -2, -1, 0, 2, 1};
1060ef8b56bfSJed Brown   static const PetscInt quadInv[8]   = {-4, -3, -2, -1, 0, 3, 2, 1};
1061ef8b56bfSJed 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};
10629371c9d4SSatish 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};
1063ef8b56bfSJed Brown   static const PetscInt tripInv[12]  = {-5, -6, -4, -3, -2, -1, 0, 2, 1, 3, 4, 5};
1064ef8b56bfSJed Brown   static const PetscInt ttriInv[12]  = {-6, -5, -4, -3, -2, -1, 0, 2, 1, 3, 5, 4};
1065ef8b56bfSJed Brown   static const PetscInt tquadInv[16] = {-8, -7, -6, -5, -4, -3, -2, -1, 0, 3, 2, 1, 4, 7, 6, 5};
1066ef8b56bfSJed Brown   static const PetscInt pyrInv[8]    = {-4, -3, -2, -1, 0, 3, 2, 1};
1067b5a892a1SMatthew G. Knepley   switch (ct) {
1068d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
1069d71ae5a4SJacob Faibussowitsch     return 0;
1070b5a892a1SMatthew G. Knepley   case DM_POLYTOPE_SEGMENT:
1071d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
1072d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, o2);
1073d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
1074d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, triInv[o2 + 3]);
1075d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
1076d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, quadInv[o2 + 4]);
1077d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
1078d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, o2);
1079d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
1080d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, tetInv[o2 + 12]);
1081d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
1082d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, hexInv[o2 + 24]);
1083d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
1084d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, tripInv[o2 + 6]);
1085d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
1086d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, ttriInv[o2 + 6]);
1087d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
1088d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, tquadInv[o2 + 8]);
1089d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
1090d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, pyrInv[o2 + 4]);
1091d71ae5a4SJacob Faibussowitsch   default:
1092d71ae5a4SJacob Faibussowitsch     return 0;
1093b5a892a1SMatthew G. Knepley   }
1094b5a892a1SMatthew G. Knepley }
1095b5a892a1SMatthew G. Knepley 
1096b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeMatchOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *);
1097b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeMatchVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *);
1098b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeGetOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *);
1099b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeGetVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *);
1100012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeInCellTest(DMPolytopeType, const PetscReal[], PetscBool *);
1101b5a892a1SMatthew G. Knepley 
110207218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTComputeFaceQuadPermutation(DMPolytopeType, PetscQuadrature, PetscInt *, IS *[]);
110307218a29SMatthew G. Knepley 
1104e1589f56SBarry Smith #endif
1105