xref: /petsc/include/petscdm.h (revision 4fb89dddf56594b92bdd2ca7e24874fafe134f45)
1e1589f56SBarry Smith /*
2e1589f56SBarry Smith       Objects to manage the interactions between the mesh data structures and the algebraic objects
3e1589f56SBarry Smith */
426bd1501SBarry Smith #if !defined(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>
11e1589f56SBarry Smith 
12ac09b921SBarry Smith /* SUBMANSEC = DM */
13ac09b921SBarry Smith 
14607a6623SBarry Smith PETSC_EXTERN PetscErrorCode DMInitializePackage(void);
15e1589f56SBarry Smith 
16014dd563SJed Brown PETSC_EXTERN PetscClassId DM_CLASSID;
17d67d17b1SMatthew G. Knepley PETSC_EXTERN PetscClassId DMLABEL_CLASSID;
18e1589f56SBarry Smith 
19528c63e0SDave May #define DMLOCATEPOINT_POINT_NOT_FOUND -367
20528c63e0SDave May 
2176bdecfbSBarry Smith /*J
228f6c3df8SBarry Smith     DMType - String with the name of a PETSc DM
23e1589f56SBarry Smith 
24e1589f56SBarry Smith    Level: beginner
25e1589f56SBarry Smith 
26db781477SPatrick Sanan .seealso: `DMSetType()`, `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[];
47140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList DMList;
48c0517cd5SMatthew G. Knepley PETSC_EXTERN DMGeneratorFunctionList DMGenerateList;
49014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreate(MPI_Comm,DM*);
5038221697SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMClone(DM,DM*);
5119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetType(DM, DMType);
5219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMGetType(DM, DMType *);
53bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode DMRegister(const char[],PetscErrorCode (*)(DM));
54014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRegisterDestroy(void);
55e1589f56SBarry Smith 
56014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMView(DM,PetscViewer);
57014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLoad(DM,PetscViewer);
58014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDestroy(DM*);
59014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*);
60014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*);
61014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM,Vec *);
62014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM,Vec *);
63014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM,Vec *);
64014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM,Vec *);
65014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
6650eeb1caSToby Isaac PETSC_EXTERN PetscErrorCode DMClearLocalVectors(DM);
67e77ac854SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasNamedGlobalVector(DM,const char*,PetscBool*);
68014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM,const char*,Vec*);
69014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM,const char*,Vec*);
70e77ac854SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasNamedLocalVector(DM,const char*,PetscBool*);
712348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM,const char*,Vec*);
722348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM,const char*,Vec*);
73014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*);
74014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM,PetscInt*,char***,IS**);
75014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM,PetscInt*);
76b412c318SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateColoring(DM,ISColoringType,ISColoring*);
77b412c318SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM,Mat*);
78aa0f6e3cSJed Brown PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateSkip(DM,PetscBool);
79014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM,PetscBool);
80b06ff27eSHong Zhang PETSC_EXTERN PetscErrorCode DMSetMatrixStructureOnly(DM,PetscBool);
81014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM,DM,Mat*,Vec*);
823ad4599aSBarry Smith PETSC_EXTERN PetscErrorCode DMCreateRestriction(DM,DM,Mat*);
836dbf9973SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMCreateInjection(DM,DM,Mat*);
84bd041c0cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateMassMatrix(DM,DM,Mat*);
85b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateMassMatrixLumped(DM,Vec*);
8669291d52SBarry Smith PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,MPI_Datatype,void*);
8769291d52SBarry Smith PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,MPI_Datatype,void*);
88014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*);
89014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*);
90a8fb8f29SToby Isaac PETSC_EXTERN PetscErrorCode DMGetCoarseDM(DM,DM*);
91a8fb8f29SToby Isaac PETSC_EXTERN PetscErrorCode DMSetCoarseDM(DM,DM);
9288bdff64SToby Isaac PETSC_EXTERN PetscErrorCode DMGetFineDM(DM,DM*);
9388bdff64SToby Isaac PETSC_EXTERN PetscErrorCode DMSetFineDM(DM,DM);
94014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]);
95014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]);
96014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
97dc822a44SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
98014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
993d8e3701SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
100014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM);
101014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM);
1021f3379b2SToby Isaac PETSC_EXTERN PetscErrorCode DMInterpolateSolution(DM,DM,Mat,Vec,Vec);
103d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMExtrude(DM,PetscInt,DM*);
104014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
105fe2efc57SMark PETSC_EXTERN PetscErrorCode DMViewFromOptions(DM,PetscObject,const char[]);
106ca266f36SBarry Smith 
107c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerate(DM, const char [], PetscBool , DM *);
1089fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMGenerateRegister(const char[],PetscErrorCode (*)(DM,PetscBool,DM*),PetscErrorCode (*)(DM,PetscReal*,DM*),PetscErrorCode (*)(DM,Vec,DMLabel,DMLabel,DM*),PetscInt);
109c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerateRegisterAll(void);
110c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerateRegisterDestroy(void);
111a1b0c543SToby Isaac PETSC_EXTERN PetscErrorCode DMAdaptLabel(DM,DMLabel,DM*);
1129fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMAdaptMetric(DM, Vec, DMLabel, DMLabel, DM *);
113df0b854cSToby Isaac 
114014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetUp(DM);
115014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM,DM,Mat,Vec*);
11697779f9aSLisandro Dalcin PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use DMDACreateAggregates() or DMCreateRestriction() (since version 3.12)") PetscErrorCode DMCreateAggregates(DM,DM,Mat*);
117baf369e7SPeter Brune PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
118d4d07f1eSToby Isaac PETSC_EXTERN PetscErrorCode DMLocalToGlobalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
11901729b5cSPatrick Sanan PETSC_EXTERN PetscErrorCode DMGlobalToLocal(DM,Vec,InsertMode,Vec);
120014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec);
121014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec);
12201729b5cSPatrick Sanan PETSC_EXTERN PetscErrorCode DMLocalToGlobal(DM,Vec,InsertMode,Vec);
123014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec);
124014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec);
125d78e899eSRichard Tran Mills PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM,Vec,InsertMode,Vec);
126d78e899eSRichard Tran Mills PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM,Vec,InsertMode,Vec);
12719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMConvert(DM,DMType,DM*);
128e1589f56SBarry Smith 
129c73cfb54SMatthew G. Knepley /* Topology support */
130c73cfb54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDimension(DM,PetscInt*);
131c73cfb54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetDimension(DM,PetscInt);
132793f3fe5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDimPoints(DM,PetscInt,PetscInt*,PetscInt*);
1338e4ac7eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetUseNatural(DM,PetscBool*);
1348e4ac7eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetUseNatural(DM,PetscBool);
1356858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNeighbors(DM,PetscInt*,const PetscMPIInt**);
13646e270d4SMatthew G. Knepley 
13746e270d4SMatthew G. Knepley /* Coordinate support */
1386636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM,DM*);
1391cfe2091SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM,DM);
1406858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinateDM(DM,DM*);
1416858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinateDM(DM,DM);
14246e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateDim(DM,PetscInt*);
14346e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateDim(DM,PetscInt);
144e8abe2deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM,PetscSection*);
14546e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM,PetscInt,PetscSection);
1466858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinateSection(DM,PetscSection*);
1476858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinateSection(DM,PetscInt,PetscSection);
1486636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM,Vec*);
1496636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM,Vec);
1506858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinates(DM,Vec*);
1516858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinates(DM,Vec);
15281e9a530SVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalSetUp(DM);
1536858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM,Vec*);
15481e9a530SVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalNoncollective(DM,Vec*);
1552db98f8dSVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalTuple(DM,IS,PetscSection*,Vec*);
1566636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM,Vec);
1576858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM);
1586858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocal(DM,Vec*);
1596858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM,Vec*);
1606858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinatesLocal(DM,Vec);
1616858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateField(DM,DMField*);
1626858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateField(DM,DMField);
1636858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetLocalBoundingBox(DM,PetscReal[],PetscReal[]);
1646858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetBoundingBox(DM,PetscReal[],PetscReal[]);
1656858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectCoordinates(DM,PetscFE);
16662a38674SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocatePoints(DM,Vec,DMPointLocationType,PetscSF*);
1676858538eSMatthew G. Knepley 
1686858538eSMatthew G. Knepley /* Periodicity support */
169*4fb89dddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM,const PetscReal*[], const PetscReal*[], const PetscReal*[]);
170*4fb89dddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM,const PetscReal[], const PetscReal[], const PetscReal[]);
171e907e85cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate(DM, const PetscScalar[], PetscBool, PetscScalar[]);
1722e17dfb7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocalizeCoordinates(DM);
17336447a5eSToby Isaac PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalized(DM,PetscBool*);
1748f700142SStefano Zampini PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalizedLocal(DM,PetscBool*);
1756636e97aSMatthew G Knepley 
1765dbd56e3SPeter Brune /* block hook interface */
177be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
178b3a6b972SJed Brown PETSC_EXTERN PetscErrorCode DMSubDomainHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
179be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM);
1805dbd56e3SPeter Brune 
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []);
18231697293SDave May PETSC_EXTERN PetscErrorCode DMAppendOptionsPrefix(DM,const char []);
18331697293SDave May PETSC_EXTERN PetscErrorCode DMGetOptionsPrefix(DM,const char*[]);
18419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType);
185c0dedaeaSBarry Smith PETSC_EXTERN PetscErrorCode DMGetVecType(DM,VecType*);
18619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType);
187c0dedaeaSBarry Smith PETSC_EXTERN PetscErrorCode DMGetMatType(DM,MatType*);
1888f1509bcSBarry Smith PETSC_EXTERN PetscErrorCode DMSetISColoringType(DM,ISColoringType);
1898f1509bcSBarry Smith PETSC_EXTERN PetscErrorCode DMGetISColoringType(DM,ISColoringType*);
190014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*);
191014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**));
192014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*);
193014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec));
194014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *);
195b0ae01b7SPeter Brune PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *);
1963ad4599aSBarry Smith PETSC_EXTERN PetscErrorCode DMHasCreateRestriction(DM,PetscBool *);
197a7058e45SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMHasCreateInjection(DM,PetscBool *);
198014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec);
19993d92d96SBarry Smith 
20037bc7515SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, const PetscInt[], IS *, DM *);
2012adcc780SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSuperDM(DM[], PetscInt, IS **, DM *);
202792b654fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSectionSubDM(DM,PetscInt,const PetscInt[],IS*,DM*);
203792b654fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSectionSuperDM(DM[],PetscInt,IS**,DM*);
20416621825SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**);
2058d4ac253SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**);
206e30e807fSPeter Brune PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
207e7c4fc90SDmitry Karpeev 
208014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*);
209fef3a512SBarry Smith PETSC_EXTERN PetscErrorCode DMSetRefineLevel(DM,PetscInt);
210014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*);
2119a64c4a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoarsenLevel(DM,PetscInt);
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);
213e1589f56SBarry Smith 
2145f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*);
2155f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
216c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*);
217c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
218531c7667SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringUseDM(Mat,MatFDColoring);
2195f1ad066SMatthew G Knepley 
220e1589f56SBarry Smith typedef struct NLF_DAAD* NLF;
221e1589f56SBarry Smith 
222bc2bf880SBarry Smith #define DM_FILE_CLASSID 1211221
2237da65231SMatthew G Knepley 
2247da65231SMatthew G Knepley /* FEM support */
225014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []);
226014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []);
2276113b454SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char [], PetscReal, Vec);
2287da65231SMatthew G Knepley 
2298cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *));
2308cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *));
2318cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *));
2328cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *));
233f9d4088aSMatthew G. Knepley 
234061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMGetSection(DM, PetscSection *); /* Use DMGetLocalSection() in new code (since v3.12) */
235061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMSetSection(DM, PetscSection);   /* Use DMSetLocalSection() in new code (since v3.12) */
236061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalSection(DM, PetscSection *);
237061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMSetLocalSection(DM, PetscSection);
238e87a4003SBarry Smith PETSC_EXTERN PetscErrorCode DMGetGlobalSection(DM, PetscSection *);
239e87a4003SBarry Smith PETSC_EXTERN PetscErrorCode DMSetGlobalSection(DM, PetscSection);
2409fbee547SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Use DMGetSection() (since v3.9)") PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *s) {return DMGetSection(dm,s);}
2419fbee547SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Use DMSetSection() (since v3.9)") PetscErrorCode DMSetDefaultSection(DM dm, PetscSection s) {return DMSetSection(dm,s);}
2429fbee547SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Use DMGetGlobalSection() (since v3.9)") PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *s) {return DMGetGlobalSection(dm,s);}
2439fbee547SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Use DMSetGlobalSection() (since v3.9)") PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection s) {return DMSetGlobalSection(dm,s);}
244e87a4003SBarry Smith 
2451bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMGetSectionSF(DM, PetscSF*);
2461bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMSetSectionSF(DM, PetscSF);
2471bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateSectionSF(DM, PetscSection, PetscSection);
2489fbee547SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Use DMGetSectionSF() (since v3.12)") PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *s) {return DMGetSectionSF(dm,s);}
2499fbee547SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Use DMSetSectionSF() (since v3.12)") PetscErrorCode DMSetDefaultSF(DM dm, PetscSF s) {return DMSetSectionSF(dm,s);}
2509fbee547SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Use DMCreateSectionSF() (since v3.12)") PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection l, PetscSection g) {return DMCreateSectionSF(dm,l,g);}
2511bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
2521bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);
2534f37162bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNaturalSF(DM, PetscSF *);
2544f37162bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNaturalSF(DM, PetscSF);
2551bb6d2a8SBarry Smith 
25679769bd5SJed Brown PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *, Vec *);
25779769bd5SJed Brown PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat, Vec);
258e87a4003SBarry Smith 
25914f150ffSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *);
260cdb7a50dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *);
261cdb7a50dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal);
262cdb7a50dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char *, PetscInt, PetscReal *);
26314f150ffSMatthew G. Knepley 
264af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
265af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
26644a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, DMLabel *, PetscObject *);
26744a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, DMLabel, PetscObject);
26844a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAddField(DM, DMLabel, PetscObject);
269e0b68406SMatthew Knepley PETSC_EXTERN PetscErrorCode DMSetFieldAvoidTensor(DM, PetscInt, PetscBool);
270e0b68406SMatthew Knepley PETSC_EXTERN PetscErrorCode DMGetFieldAvoidTensor(DM, PetscInt, PetscBool *);
27144a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMClearFields(DM);
272e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyFields(DM, DM);
27334aa8a36SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAdjacency(DM, PetscInt, PetscBool *, PetscBool *);
27434aa8a36SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetAdjacency(DM, PetscInt, PetscBool, PetscBool);
275b0441da4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetBasicAdjacency(DM, PetscBool *, PetscBool *);
276b0441da4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetBasicAdjacency(DM, PetscBool, PetscBool);
277e5e52638SMatthew G. Knepley 
278e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNumDS(DM, PetscInt *);
279e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *);
280e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellDS(DM, PetscInt, PetscDS *);
281b3cf3223SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetRegionDS(DM, DMLabel, IS *, PetscDS *);
282b3cf3223SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetRegionDS(DM, DMLabel, IS, PetscDS);
283b3cf3223SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetRegionNumDS(DM, PetscInt, DMLabel *, IS *, PetscDS *);
284083401c6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetRegionNumDS(DM, PetscInt, DMLabel, IS, PetscDS);
2851d3af9e0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMFindRegionNum(DM, PetscDS, PetscInt *);
2862df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateFEDefault(DM, PetscInt, const char[], PetscInt, PetscFE *);
287e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateDS(DM);
288e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMClearDS(DM);
289e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyDS(DM, DM);
290e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyDisc(DM, DM);
291f2cacb80SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMComputeExactSolution(DM, PetscReal, Vec, Vec);
2929a2a23afSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNumAuxiliaryVec(DM, PetscInt *);
293ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAuxiliaryVec(DM, DMLabel, PetscInt, PetscInt, Vec *);
294ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetAuxiliaryVec(DM, DMLabel, PetscInt, PetscInt, Vec);
295ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAuxiliaryLabels(DM, DMLabel[], PetscInt[], PetscInt[]);
2969a2a23afSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyAuxiliaryVec(DM, DM);
297af122d2aSMatthew G Knepley 
2984267b1a3SMatthew G. Knepley /*MC
2994267b1a3SMatthew G. Knepley   DMInterpolationInfo - Structure for holding information about interpolation on a mesh
3004267b1a3SMatthew G. Knepley 
3014267b1a3SMatthew G. Knepley   Level: intermediate
3024267b1a3SMatthew G. Knepley 
3034267b1a3SMatthew G. Knepley   Synopsis:
3044267b1a3SMatthew G. Knepley     comm   - The communicator
3054267b1a3SMatthew G. Knepley     dim    - The spatial dimension of points
3064267b1a3SMatthew G. Knepley     nInput - The number of input points
3074267b1a3SMatthew G. Knepley     points - The input point coordinates
3084267b1a3SMatthew G. Knepley     cells  - The cell containing each point
3094267b1a3SMatthew G. Knepley     n      - The number of local points
3104267b1a3SMatthew G. Knepley     coords - The point coordinates
3114267b1a3SMatthew G. Knepley     dof    - The number of components to interpolate
3124267b1a3SMatthew G. Knepley 
313db781477SPatrick Sanan .seealso: `DMInterpolationCreate()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
3144267b1a3SMatthew G. Knepley M*/
315e87bb0d3SMatthew G Knepley struct _DMInterpolationInfo {
316e87bb0d3SMatthew G Knepley   MPI_Comm   comm;
317e87bb0d3SMatthew G Knepley   PetscInt   dim;    /*1 The spatial dimension of points */
318e87bb0d3SMatthew G Knepley   PetscInt   nInput; /* The number of input points */
319e87bb0d3SMatthew G Knepley   PetscReal *points; /* The input point coordinates */
320e87bb0d3SMatthew G Knepley   PetscInt  *cells;  /* The cell containing each point */
321e87bb0d3SMatthew G Knepley   PetscInt   n;      /* The number of local points */
322e87bb0d3SMatthew G Knepley   Vec        coords; /* The point coordinates */
323e87bb0d3SMatthew G Knepley   PetscInt   dof;    /* The number of components to interpolate */
324e87bb0d3SMatthew G Knepley };
325e87bb0d3SMatthew G Knepley typedef struct _DMInterpolationInfo *DMInterpolationInfo;
326e87bb0d3SMatthew G Knepley 
32794b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
32894b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
32994b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
33094b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
33194b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
33294b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
33352aa1562SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool, PetscBool);
33494b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
33594b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
33694b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
33794b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
33894b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);
339c58f1c22SToby Isaac 
340c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMCreateLabel(DM, const char []);
341c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelValue(DM, const char[], PetscInt, PetscInt *);
342c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMSetLabelValue(DM, const char[], PetscInt, PetscInt);
343c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMClearLabelValue(DM, const char[], PetscInt, PetscInt);
344c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelSize(DM, const char[], PetscInt *);
345c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelIdIS(DM, const char[], IS *);
346c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetStratumSize(DM, const char [], PetscInt, PetscInt *);
347c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetStratumIS(DM, const char [], PetscInt, IS *);
3484de306b1SToby Isaac PETSC_EXTERN PetscErrorCode DMSetStratumIS(DM, const char [], PetscInt, IS);
349c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMClearLabelStratum(DM, const char[], PetscInt);
350c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelOutput(DM, const char[], PetscBool *);
351c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMSetLabelOutput(DM, const char[], PetscBool);
3525f15299fSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMGetFirstLabeledPoint(DM, DM, DMLabel, PetscInt, const PetscInt *, PetscInt, PetscInt *, PetscDS *);
353c58f1c22SToby Isaac 
3542cbb9b06SVaclav Hapla /*E
3552cbb9b06SVaclav Hapla    DMCopyLabelsMode - Determines how DMCopyLabels() behaves when there is a DMLabel in the source and destination DMs with the same name
3562cbb9b06SVaclav Hapla 
3572cbb9b06SVaclav Hapla    Level: advanced
3582cbb9b06SVaclav Hapla 
3592cbb9b06SVaclav Hapla $ DM_COPY_LABELS_REPLACE  - replace label in destination by label from source
3602cbb9b06SVaclav Hapla $ DM_COPY_LABELS_KEEP     - keep destination label
3612cbb9b06SVaclav Hapla $ DM_COPY_LABELS_FAIL     - throw error
3622cbb9b06SVaclav Hapla 
3632cbb9b06SVaclav Hapla E*/
3642cbb9b06SVaclav Hapla typedef enum {DM_COPY_LABELS_REPLACE, DM_COPY_LABELS_KEEP, DM_COPY_LABELS_FAIL} DMCopyLabelsMode;
3652cbb9b06SVaclav Hapla PETSC_EXTERN const char *const DMCopyLabelsModes[];
3662cbb9b06SVaclav Hapla 
367c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetNumLabels(DM, PetscInt *);
368c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelName(DM, PetscInt, const char **);
369c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMHasLabel(DM, const char [], PetscBool *);
370c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabel(DM, const char *, DMLabel *);
3714a7ee7d0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetLabel(DM, DMLabel);
372c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelByNum(DM, PetscInt, DMLabel *);
373c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMAddLabel(DM, DMLabel);
374c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMRemoveLabel(DM, const char [], DMLabel *);
375306894acSVaclav Hapla PETSC_EXTERN PetscErrorCode DMRemoveLabelBySelf(DM, DMLabel *, PetscBool);
3762cbb9b06SVaclav Hapla PETSC_EXTERN PetscErrorCode DMCopyLabels(DM, DM, PetscCopyMode, PetscBool, DMCopyLabelsMode emode);
377609dae6eSVaclav Hapla PETSC_EXTERN PetscErrorCode DMCompareLabels(DM, DM, PetscBool *, char **);
378c58f1c22SToby Isaac 
37945480ffeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAddBoundary(DM, DMBoundaryConditionType, const char[], DMLabel, PetscInt, const PetscInt[], PetscInt, PetscInt, const PetscInt[], void (*)(void), void (*)(void), void *, PetscInt *);
380a6ba4734SToby Isaac PETSC_EXTERN PetscErrorCode DMIsBoundaryPoint(DM, PetscInt, PetscBool *);
3814d6f44ffSToby Isaac 
3820709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMProjectFunction(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
3830709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMProjectFunctionLocal(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
3842c53366bSMatthew 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);
3851c531cf8SMatthew 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);
386191494d9SMatthew 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);
3871c531cf8SMatthew 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);
388ece3a9fcSMatthew 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);
3890709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2Diff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
390b698f381SToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2GradientDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal [], PetscReal *);
3911189c1efSToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2FieldDiff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
3922e4af2aeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMComputeError(DM, Vec, PetscReal[], Vec *);
393ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasBasisTransform(DM,PetscBool*);
394ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyTransform(DM, DM);
3958320bc6fSPatrick Sanan 
3968320bc6fSPatrick Sanan PETSC_EXTERN PetscErrorCode DMGetCompatibility(DM,DM,PetscBool*,PetscBool*);
397c0f0dcc3SMatthew G. Knepley 
398c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitorSet(DM, PetscErrorCode (*)(DM, void *), void *, PetscErrorCode (*)(void**));
399c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitorCancel(DM);
400c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitorSetFromOptions(DM, const char[], const char[], const char[], PetscErrorCode (*)(DM, void *), PetscErrorCode (*)(DM, PetscViewerAndFormat *), PetscBool *);
401c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitor(DM);
402c0f0dcc3SMatthew G. Knepley 
4039fbee547SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetDim(DMPolytopeType ct)
40470a7d78aSStefano Zampini {
405412e9a14SMatthew G. Knepley   switch (ct) {
406412e9a14SMatthew G. Knepley     case DM_POLYTOPE_POINT:
407412e9a14SMatthew G. Knepley       return 0;
408412e9a14SMatthew G. Knepley     case DM_POLYTOPE_SEGMENT:
409412e9a14SMatthew G. Knepley     case DM_POLYTOPE_POINT_PRISM_TENSOR:
410412e9a14SMatthew G. Knepley       return 1;
411412e9a14SMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:
412412e9a14SMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:
413412e9a14SMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR:
414412e9a14SMatthew G. Knepley       return 2;
415412e9a14SMatthew G. Knepley     case DM_POLYTOPE_TETRAHEDRON:
416412e9a14SMatthew G. Knepley     case DM_POLYTOPE_HEXAHEDRON:
417412e9a14SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM:
418412e9a14SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM_TENSOR:
419412e9a14SMatthew G. Knepley     case DM_POLYTOPE_QUAD_PRISM_TENSOR:
420da9060c4SMatthew G. Knepley     case DM_POLYTOPE_PYRAMID:
421412e9a14SMatthew G. Knepley       return 3;
422412e9a14SMatthew G. Knepley     default: return -1;
423412e9a14SMatthew G. Knepley   }
424412e9a14SMatthew G. Knepley }
425412e9a14SMatthew G. Knepley 
4269fbee547SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetConeSize(DMPolytopeType ct)
427412e9a14SMatthew G. Knepley {
428412e9a14SMatthew G. Knepley   switch (ct) {
429412e9a14SMatthew G. Knepley     case DM_POLYTOPE_POINT:              return 0;
430412e9a14SMatthew G. Knepley     case DM_POLYTOPE_SEGMENT:            return 2;
431412e9a14SMatthew G. Knepley     case DM_POLYTOPE_POINT_PRISM_TENSOR: return 2;
432412e9a14SMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:           return 3;
433412e9a14SMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:      return 4;
434412e9a14SMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return 4;
435412e9a14SMatthew G. Knepley     case DM_POLYTOPE_TETRAHEDRON:        return 4;
436412e9a14SMatthew G. Knepley     case DM_POLYTOPE_HEXAHEDRON:         return 6;
437412e9a14SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM:          return 5;
438412e9a14SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return 5;
439412e9a14SMatthew G. Knepley     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return 6;
440da9060c4SMatthew G. Knepley     case DM_POLYTOPE_PYRAMID:            return 5;
441412e9a14SMatthew G. Knepley     default: return -1;
442412e9a14SMatthew G. Knepley   }
443412e9a14SMatthew G. Knepley }
444412e9a14SMatthew G. Knepley 
4459fbee547SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetNumVertices(DMPolytopeType ct)
446412e9a14SMatthew G. Knepley {
447412e9a14SMatthew G. Knepley   switch (ct) {
448412e9a14SMatthew G. Knepley     case DM_POLYTOPE_POINT:              return 1;
449412e9a14SMatthew G. Knepley     case DM_POLYTOPE_SEGMENT:            return 2;
450412e9a14SMatthew G. Knepley     case DM_POLYTOPE_POINT_PRISM_TENSOR: return 2;
451412e9a14SMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:           return 3;
452412e9a14SMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:      return 4;
453412e9a14SMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return 4;
454412e9a14SMatthew G. Knepley     case DM_POLYTOPE_TETRAHEDRON:        return 4;
455412e9a14SMatthew G. Knepley     case DM_POLYTOPE_HEXAHEDRON:         return 8;
456412e9a14SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM:          return 6;
457412e9a14SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return 6;
458412e9a14SMatthew G. Knepley     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return 8;
459da9060c4SMatthew G. Knepley     case DM_POLYTOPE_PYRAMID:            return 5;
460412e9a14SMatthew G. Knepley     default: return -1;
461412e9a14SMatthew G. Knepley   }
462412e9a14SMatthew G. Knepley }
463412e9a14SMatthew G. Knepley 
4649fbee547SJacob Faibussowitsch static inline DMPolytopeType DMPolytopeTypeSimpleShape(PetscInt dim, PetscBool simplex)
4659318fe57SMatthew G. Knepley {
4669318fe57SMatthew G. Knepley   return dim == 0 ? DM_POLYTOPE_POINT :
4679318fe57SMatthew G. Knepley         (dim == 1 ? DM_POLYTOPE_SEGMENT :
4689318fe57SMatthew G. Knepley         (dim == 2 ? (simplex ? DM_POLYTOPE_TRIANGLE : DM_POLYTOPE_QUADRILATERAL) :
4699318fe57SMatthew G. Knepley         (dim == 3 ? (simplex ? DM_POLYTOPE_TETRAHEDRON : DM_POLYTOPE_HEXAHEDRON) : DM_POLYTOPE_UNKNOWN)));
4709318fe57SMatthew G. Knepley }
4719318fe57SMatthew G. Knepley 
4729fbee547SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetNumArrangments(DMPolytopeType ct)
473b5a892a1SMatthew G. Knepley {
474b5a892a1SMatthew G. Knepley   switch (ct) {
475b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_POINT:              return 1;
476b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_SEGMENT:            return 2;
477b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_POINT_PRISM_TENSOR: return 2;
478b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:           return 6;
479b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:      return 8;
480b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return 4;
481b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TETRAHEDRON:        return 24;
482b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_HEXAHEDRON:         return 48;
483b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM:          return 12;
484b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return 12;
485b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return 16;
486b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_PYRAMID:            return 8;
487b5a892a1SMatthew G. Knepley     default: return -1;
488b5a892a1SMatthew G. Knepley   }
489b5a892a1SMatthew G. Knepley }
490b5a892a1SMatthew G. Knepley 
491b5a892a1SMatthew G. Knepley /* An arrangement is a face order combined with an orientation for each face */
4929fbee547SJacob Faibussowitsch static inline const PetscInt *DMPolytopeTypeGetArrangment(DMPolytopeType ct, PetscInt o)
493b5a892a1SMatthew G. Knepley {
494ef8b56bfSJed Brown   static const PetscInt pntArr[1*2] = {0, 0};
495b5a892a1SMatthew G. Knepley   /* a: swap */
496ef8b56bfSJed Brown   static const PetscInt segArr[2*2*2] = {
497b5a892a1SMatthew G. Knepley     1, 0,  0, 0, /* -1: a */
498b5a892a1SMatthew G. Knepley     0, 0,  1, 0, /*  0: e */};
499b5a892a1SMatthew G. Knepley   /* a: swap first two
500b5a892a1SMatthew G. Knepley      b: swap last two */
501ef8b56bfSJed Brown   static const PetscInt triArr[6*3*2] = {
502b5a892a1SMatthew G. Knepley     0, -1,  2, -1,  1, -1, /* -3: b */
503b5a892a1SMatthew G. Knepley     2, -1,  1, -1,  0, -1, /* -2: aba */
504b5a892a1SMatthew G. Knepley     1, -1,  0, -1,  2, -1, /* -1: a */
505b5a892a1SMatthew G. Knepley     0,  0,  1,  0,  2,  0, /*  0: identity */
506b5a892a1SMatthew G. Knepley     1,  0,  2,  0,  0,  0, /*  1: ba */
507b5a892a1SMatthew G. Knepley     2,  0,  0,  0,  1,  0, /*  2: ab */};
508b5a892a1SMatthew G. Knepley   /* a: forward cyclic permutation
509b5a892a1SMatthew G. Knepley      b: swap first and last pairs */
510ef8b56bfSJed Brown   static const PetscInt quadArr[8*4*2] = {
511b5a892a1SMatthew G. Knepley     1, -1,  0, -1,  3, -1,  2, -1, /* -4: b */
512b5a892a1SMatthew G. Knepley     0, -1,  3, -1,  2, -1,  1, -1, /* -3: b a^3 = a b */
513b5a892a1SMatthew G. Knepley     3, -1,  2, -1,  1, -1,  0, -1, /* -2: b a^2 = a^2 b */
514b5a892a1SMatthew G. Knepley     2, -1,  1, -1,  0, -1,  3, -1, /* -1: b a   = a^3 b */
515b5a892a1SMatthew G. Knepley     0,  0,  1,  0,  2,  0,  3,  0, /*  0: identity */
516b5a892a1SMatthew G. Knepley     1,  0,  2,  0,  3,  0,  0,  0, /*  1: a */
517b5a892a1SMatthew G. Knepley     2,  0,  3,  0,  0,  0,  1,  0, /*  2: a^2 */
518b5a892a1SMatthew G. Knepley     3,  0,  0,  0,  1,  0,  2,  0, /*  3: a^3 */};
519b5a892a1SMatthew G. Knepley   /* r: rotate 180
520b5a892a1SMatthew G. Knepley      b: swap top and bottom segments */
521ef8b56bfSJed Brown   static const PetscInt tsegArr[4*4*2] = {
522b5a892a1SMatthew G. Knepley     1, -1,  0, -1,  3, -1,  2, -1, /* -2: r b */
523b5a892a1SMatthew G. Knepley     0, -1,  1, -1,  3,  0,  2,  0, /* -1: r */
524b5a892a1SMatthew G. Knepley     0,  0,  1,  0,  2,  0,  3,  0, /*  0: identity */
525b5a892a1SMatthew G. Knepley     1,  0,  0,  0,  2, -1,  3, -1, /*  1: b */};
526b5a892a1SMatthew G. Knepley   /* https://en.wikiversity.org/wiki/Symmetric_group_S4 */
527ef8b56bfSJed Brown   static const PetscInt tetArr[24*4*2] = {
528b5a892a1SMatthew G. Knepley     3, -2,  2, -3,  0, -1,  1, -1, /* -12: (1324)   p22 */
529b5a892a1SMatthew G. Knepley     3, -1,  1, -3,  2, -1,  0, -1, /* -11: (14)     p21 */
530b5a892a1SMatthew G. Knepley     3, -3,  0, -3,  1, -1,  2, -1, /* -10: (1234)   p18 */
531b5a892a1SMatthew G. Knepley     2, -1,  3, -1,  1, -3,  0, -2, /*  -9: (1423)   p17 */
532b5a892a1SMatthew G. Knepley     2, -3,  0, -1,  3, -2,  1, -3, /*  -8: (1342)   p13 */
533b5a892a1SMatthew G. Knepley     2, -2,  1, -2,  0, -2,  3, -2, /*  -7: (24)     p14 */
534b5a892a1SMatthew G. Knepley     1, -2,  0, -2,  2, -2,  3, -1, /*  -6: (34)     p6  */
535b5a892a1SMatthew G. Knepley     1, -1,  3, -3,  0, -3,  2, -2, /*  -5: (1243)   p10 */
536b5a892a1SMatthew G. Knepley     1, -3,  2, -1,  3, -1,  0, -3, /*  -4: (1432)   p9  */
537b5a892a1SMatthew G. Knepley     0, -3,  1, -1,  3, -3,  2, -3, /*  -3: (12)     p1  */
538b5a892a1SMatthew G. Knepley     0, -2,  2, -2,  1, -2,  3, -3, /*  -2: (23)     p2  */
539b5a892a1SMatthew G. Knepley     0, -1,  3, -2,  2, -3,  1, -2, /*  -1: (13)     p5  */
540b5a892a1SMatthew G. Knepley     0,  0,  1,  0,  2,  0,  3,  0, /*   0: ()       p0  */
541b5a892a1SMatthew G. Knepley     0,  1,  3,  1,  1,  2,  2,  0, /*   1: (123)    p4  */
542b5a892a1SMatthew G. Knepley     0,  2,  2,  1,  3,  0,  1,  2, /*   2: (132)    p3  */
543b5a892a1SMatthew G. Knepley     1,  2,  0,  1,  3,  1,  2,  2, /*   3: (12)(34) p7  */
544b5a892a1SMatthew G. Knepley     1,  0,  2,  0,  0,  0,  3,  1, /*   4: (243)    p8  */
545b5a892a1SMatthew G. Knepley     1,  1,  3,  2,  2,  2,  0,  0, /*   5: (143)    p11 */
546b5a892a1SMatthew G. Knepley     2,  1,  3,  0,  0,  2,  1,  0, /*   6: (13)(24) p16 */
547b5a892a1SMatthew G. Knepley     2,  2,  1,  1,  3,  2,  0,  2, /*   7: (142)    p15 */
548b5a892a1SMatthew G. Knepley     2,  0,  0,  0,  1,  0,  3,  2, /*   8: (234)    p12 */
549b5a892a1SMatthew G. Knepley     3,  2,  2,  2,  1,  1,  0,  1, /*   9: (14)(23) p23 */
550b5a892a1SMatthew G. Knepley     3,  0,  0,  2,  2,  1,  1,  1, /*  10: (134)    p19 */
551b5a892a1SMatthew G. Knepley     3,  1,  1,  2,  0,  1,  2,  1  /*  11: (124)    p20 */};
552b5a892a1SMatthew G. Knepley   /* Each rotation determines a permutation of the four diagonals, and this defines the isomorphism with S_4 */
553ef8b56bfSJed Brown   static const PetscInt hexArr[48*6*2] = {
554b5a892a1SMatthew G. Knepley     2, -3,  3, -2,  4, -2,  5, -3,  1, -3,  0, -1, /* -24: reflect bottom and use -3 on top */
555b5a892a1SMatthew G. Knepley     4, -2,  5, -2,  0, -1,  1, -4,  3, -2,  2, -3, /* -23: reflect bottom and use -3 on top */
556b5a892a1SMatthew G. Knepley     5, -3,  4, -1,  1, -2,  0, -3,  3, -4,  2, -1, /* -22: reflect bottom and use -3 on top */
557b5a892a1SMatthew G. Knepley     3, -1,  2, -4,  4, -4,  5, -1,  0, -4,  1, -4, /* -21: reflect bottom and use -3 on top */
558b5a892a1SMatthew G. Knepley     3, -3,  2, -2,  5, -1,  4, -4,  1, -1,  0, -3, /* -20: reflect bottom and use -3 on top */
559b5a892a1SMatthew G. Knepley     4, -4,  5, -4,  1, -4,  0, -1,  2, -4,  3, -1, /* -19: reflect bottom and use -3 on top */
560b5a892a1SMatthew G. Knepley     2, -1,  3, -4,  5, -3,  4, -2,  0, -2,  1, -2, /* -18: reflect bottom and use -3 on top */
561b5a892a1SMatthew G. Knepley     5, -1,  4, -3,  0, -3,  1, -2,  2, -2,  3, -3, /* -17: reflect bottom and use -3 on top */
562b5a892a1SMatthew G. Knepley     4, -3,  5, -1,  3, -2,  2, -4,  1, -4,  0, -4, /* -16: reflect bottom and use -3 on top */
563b5a892a1SMatthew G. Knepley     5, -4,  4, -4,  3, -4,  2, -2,  0, -3,  1, -1, /* -15: reflect bottom and use -3 on top */
564b5a892a1SMatthew G. Knepley     3, -4,  2, -1,  1, -1,  0, -4,  4, -4,  5, -4, /* -14: reflect bottom and use -3 on top */
565b5a892a1SMatthew G. Knepley     2, -2,  3, -3,  0, -2,  1, -3,  4, -2,  5, -2, /* -13: reflect bottom and use -3 on top */
566b5a892a1SMatthew G. Knepley     1, -3,  0, -1,  4, -1,  5, -4,  3, -1,  2, -4, /* -12: reflect bottom and use -3 on top */
567b5a892a1SMatthew G. Knepley     1, -1,  0, -3,  5, -4,  4, -1,  2, -1,  3, -4, /* -11: reflect bottom and use -3 on top */
568b5a892a1SMatthew G. Knepley     5, -2,  4, -2,  2, -2,  3, -4,  1, -2,  0, -2, /* -10: reflect bottom and use -3 on top */
569b5a892a1SMatthew G. Knepley     1, -2,  0, -2,  2, -1,  3, -1,  4, -1,  5, -3, /*  -9: reflect bottom and use -3 on top */
570b5a892a1SMatthew G. Knepley     4, -1,  5, -3,  2, -4,  3, -2,  0, -1,  1, -3, /*  -8: reflect bottom and use -3 on top */
571b5a892a1SMatthew G. Knepley     3, -2,  2, -3,  0, -4,  1, -1,  5, -1,  4, -3, /*  -7: reflect bottom and use -3 on top */
572b5a892a1SMatthew G. Knepley     1, -4,  0, -4,  3, -1,  2, -1,  5, -4,  4, -4, /*  -6: reflect bottom and use -3 on top */
573b5a892a1SMatthew G. Knepley     2, -4,  3, -1,  1, -3,  0, -2,  5, -3,  4, -1, /*  -5: reflect bottom and use -3 on top */
574b5a892a1SMatthew G. Knepley     0, -4,  1, -4,  4, -3,  5, -2,  2, -3,  3, -2, /*  -4: reflect bottom and use -3 on top */
575b5a892a1SMatthew G. Knepley     0, -3,  1, -1,  3, -3,  2, -3,  4, -3,  5, -1, /*  -3: reflect bottom and use -3 on top */
576b5a892a1SMatthew G. Knepley     0, -2,  1, -2,  5, -2,  4, -3,  3, -3,  2, -2, /*  -2: reflect bottom and use -3 on top */
577b5a892a1SMatthew G. Knepley     0, -1,  1, -3,  2, -3,  3, -3,  5, -2,  4, -2, /*  -1: reflect bottom and use -3 on top */
578b5a892a1SMatthew G. Knepley     0,  0,  1,  0,  2,  0,  3,  0,  4,  0,  5,  0, /*   0: identity */
579b5a892a1SMatthew G. Knepley     0,  1,  1,  3,  5,  3,  4,  0,  2,  0,  3,  1, /*   1: 90  rotation about z */
580b5a892a1SMatthew G. Knepley     0,  2,  1,  2,  3,  0,  2,  0,  5,  3,  4,  1, /*   2: 180 rotation about z */
581b5a892a1SMatthew G. Knepley     0,  3,  1,  1,  4,  0,  5,  3,  3,  0,  2,  1, /*   3: 270 rotation about z */
582b5a892a1SMatthew G. Knepley     2,  3,  3,  2,  1,  0,  0,  3,  4,  3,  5,  1, /*   4: 90  rotation about x */
583b5a892a1SMatthew G. Knepley     1,  3,  0,  1,  3,  2,  2,  2,  4,  2,  5,  2, /*   5: 180 rotation about x */
584b5a892a1SMatthew G. Knepley     3,  1,  2,  0,  0,  1,  1,  2,  4,  1,  5,  3, /*   6: 270 rotation about x */
585b5a892a1SMatthew G. Knepley     4,  0,  5,  0,  2,  1,  3,  3,  1,  1,  0,  3, /*   7: 90  rotation about y */
586b5a892a1SMatthew G. Knepley     1,  1,  0,  3,  2,  2,  3,  2,  5,  1,  4,  3, /*   8: 180 rotation about y */
587b5a892a1SMatthew G. Knepley     5,  1,  4,  3,  2,  3,  3,  1,  0,  0,  1,  0, /*   9: 270 rotation about y */
588b5a892a1SMatthew G. Knepley     1,  0,  0,  0,  5,  1,  4,  2,  3,  2,  2,  3, /*  10: 180 rotation about x+y */
589b5a892a1SMatthew G. Knepley     1,  2,  0,  2,  4,  2,  5,  1,  2,  2,  3,  3, /*  11: 180 rotation about x-y */
590b5a892a1SMatthew G. Knepley     2,  1,  3,  0,  0,  3,  1,  0,  5,  0,  4,  0, /*  12: 180 rotation about y+z */
591b5a892a1SMatthew G. Knepley     3,  3,  2,  2,  1,  2,  0,  1,  5,  2,  4,  2, /*  13: 180 rotation about y-z */
592b5a892a1SMatthew G. Knepley     5,  3,  4,  1,  3,  1,  2,  3,  1,  3,  0,  1, /*  14: 180 rotation about z+x */
593b5a892a1SMatthew G. Knepley     4,  2,  5,  2,  3,  3,  2,  1,  0,  2,  1,  2, /*  15: 180 rotation about z-x */
594b5a892a1SMatthew G. Knepley     5,  0,  4,  0,  0,  0,  1,  3,  3,  1,  2,  0, /*  16: 120 rotation about x+y+z (v0v6) */
595b5a892a1SMatthew G. Knepley     2,  0,  3,  1,  5,  0,  4,  3,  1,  0,  0,  0, /*  17: 240 rotation about x+y+z (v0v6) */
596b5a892a1SMatthew G. Knepley     4,  3,  5,  1,  1,  1,  0,  2,  3,  3,  2,  2, /*  18: 120 rotation about x+y-z (v4v2) */
597b5a892a1SMatthew G. Knepley     3,  2,  2,  3,  5,  2,  4,  1,  0,  1,  1,  3, /*  19: 240 rotation about x+y-z (v4v2) */
598b5a892a1SMatthew G. Knepley     3,  0,  2,  1,  4,  1,  5,  2,  1,  2,  0,  2, /*  20: 120 rotation about x-y+z (v1v5) */
599b5a892a1SMatthew G. Knepley     5,  2,  4,  2,  1,  3,  0,  0,  2,  3,  3,  2, /*  21: 240 rotation about x-y+z (v1v5) */
600b5a892a1SMatthew G. Knepley     4,  1,  5,  3,  0,  2,  1,  1,  2,  1,  3,  0, /*  22: 120 rotation about x-y-z (v7v3) */
601b5a892a1SMatthew G. Knepley     2,  2,  3,  3,  4,  3,  5,  0,  0,  3,  1,  1, /*  23: 240 rotation about x-y-z (v7v3) */
602b5a892a1SMatthew G. Knepley   };
603ef8b56bfSJed Brown   static const PetscInt tripArr[12*5*2] = {
604b5a892a1SMatthew G. Knepley     1, -3,  0, -1,  3, -1,  4, -1,  2, -1, /* -6: reflect bottom and top */
605b5a892a1SMatthew G. Knepley     1, -1,  0, -3,  4, -1,  2, -1,  3, -1, /* -5: reflect bottom and top */
606b5a892a1SMatthew G. Knepley     1, -2,  0, -2,  2, -1,  3, -1,  4, -1, /* -4: reflect bottom and top */
607b5a892a1SMatthew G. Knepley     0, -3,  1, -1,  3, -3,  2, -3,  4, -3, /* -3: reflect bottom and top */
608b5a892a1SMatthew G. Knepley     0, -2,  1, -2,  4, -3,  3, -3,  2, -3, /* -2: reflect bottom and top */
609b5a892a1SMatthew G. Knepley     0, -1,  1, -3,  2, -3,  4, -3,  3, -3, /* -1: reflect bottom and top */
610b5a892a1SMatthew G. Knepley     0,  0,  1,  0,  2,  0,  3,  0,  4,  0, /*  0: identity */
611b5a892a1SMatthew G. Knepley     0,  1,  1,  2,  4,  0,  2,  0,  3,  0, /*  1: 120 rotation about z */
612b5a892a1SMatthew G. Knepley     0,  2,  1,  1,  3,  0,  4,  0,  2,  0, /*  2: 240 rotation about z */
613b5a892a1SMatthew G. Knepley     1,  1,  0,  2,  2,  2,  4,  2,  3,  2, /*  3: 180 rotation about y of 0 */
614b5a892a1SMatthew G. Knepley     1,  0,  0,  0,  4,  2,  3,  2,  2,  2, /*  4: 180 rotation about y of 1 */
615b5a892a1SMatthew G. Knepley     1,  2,  0,  1,  3,  2,  2,  2,  4,  2, /*  5: 180 rotation about y of 2 */
616b5a892a1SMatthew G. Knepley   };
617b5a892a1SMatthew G. Knepley   /* a: rotate 120 about z
618b5a892a1SMatthew G. Knepley      b: swap top and bottom segments
619b5a892a1SMatthew G. Knepley      r: reflect */
620ef8b56bfSJed Brown   static const PetscInt ttriArr[12*5*2] = {
621b5a892a1SMatthew G. Knepley     1, -3,  0, -3,  2, -2,  4, -2,  3, -2, /* -6: r b a^2 */
622b5a892a1SMatthew G. Knepley     1, -2,  0, -2,  4, -2,  3, -2,  2, -2, /* -5: r b a */
623b5a892a1SMatthew G. Knepley     1, -1,  0, -1,  3, -2,  2, -2,  4, -2, /* -4: r b */
624b5a892a1SMatthew G. Knepley     0, -3,  1, -3,  2, -1,  4, -1,  3, -1, /* -3: r a^2 */
625b5a892a1SMatthew G. Knepley     0, -2,  1, -2,  4, -1,  3, -1,  2, -1, /* -2: r a */
626b5a892a1SMatthew G. Knepley     0, -1,  1, -1,  3, -1,  2, -1,  4, -1, /* -1: r */
627b5a892a1SMatthew G. Knepley     0,  0,  1,  0,  2,  0,  3,  0,  4,  0, /*  0: identity */
628b5a892a1SMatthew G. Knepley     0,  1,  1,  1,  3,  0,  4,  0,  2,  0, /*  1: a */
629b5a892a1SMatthew G. Knepley     0,  2,  1,  2,  4,  0,  2,  0,  3,  0, /*  2: a^2 */
630b5a892a1SMatthew G. Knepley     1,  0,  0,  0,  2,  1,  3,  1,  4,  1, /*  3: b */
631b5a892a1SMatthew G. Knepley     1,  1,  0,  1,  3,  1,  4,  1,  2,  1, /*  4: b a */
632b5a892a1SMatthew G. Knepley     1,  2,  0,  2,  4,  1,  2,  1,  3,  1, /*  5: b a^2 */
633b5a892a1SMatthew G. Knepley   };
634b5a892a1SMatthew G. Knepley   /* a: rotate 90 about z
635b5a892a1SMatthew G. Knepley      b: swap top and bottom segments
636b5a892a1SMatthew G. Knepley      r: reflect */
637ef8b56bfSJed Brown   static const PetscInt tquadArr[16*6*2] = {
638b5a892a1SMatthew G. Knepley     1, -4,  0, -4,  3, -2,  2, -2,  5, -2,  4, -2, /* -8: r b a^3 */
639b5a892a1SMatthew G. Knepley     1, -3,  0, -3,  2, -2,  5, -2,  4, -2,  3, -2, /* -7: r b a^2 */
640b5a892a1SMatthew G. Knepley     1, -2,  0, -2,  5, -2,  4, -2,  3, -2,  2, -2, /* -6: r b a */
641b5a892a1SMatthew G. Knepley     1, -1,  0, -1,  4, -2,  3, -2,  2, -2,  5, -2, /* -5: r b */
642b5a892a1SMatthew G. Knepley     0, -4,  1, -4,  3, -1,  2, -1,  5, -1,  4, -1, /* -4: r a^3 */
643b5a892a1SMatthew G. Knepley     0, -3,  1, -3,  2, -1,  5, -1,  4, -1,  3, -1, /* -3: r a^2 */
644b5a892a1SMatthew G. Knepley     0, -2,  1, -2,  5, -1,  4, -1,  3, -1,  2, -1, /* -2: r a */
645b5a892a1SMatthew G. Knepley     0, -1,  1, -1,  4, -1,  3, -1,  2, -1,  5, -1, /* -1: r */
646b5a892a1SMatthew G. Knepley     0,  0,  1,  0,  2,  0,  3,  0,  4,  0,  5,  0, /*  0: identity */
647b5a892a1SMatthew G. Knepley     0,  1,  1,  1,  3,  0,  4,  0,  5,  0,  2,  0, /*  1: a */
648b5a892a1SMatthew G. Knepley     0,  2,  1,  2,  4,  0,  5,  0,  2,  0,  3,  0, /*  2: a^2 */
649b5a892a1SMatthew G. Knepley     0,  3,  1,  3,  5,  0,  2,  0,  3,  0,  4,  0, /*  3: a^3 */
650b5a892a1SMatthew G. Knepley     1,  0,  0,  0,  2,  1,  3,  1,  4,  1,  5,  1, /*  4: b */
651b5a892a1SMatthew G. Knepley     1,  1,  0,  1,  3,  1,  4,  1,  5,  1,  2,  1, /*  5: b a */
652b5a892a1SMatthew G. Knepley     1,  2,  0,  2,  4,  1,  5,  1,  2,  1,  3,  1, /*  6: b a^2 */
653b5a892a1SMatthew G. Knepley     1,  3,  0,  3,  5,  1,  2,  1,  3,  1,  4,  1, /*  7: b a^3 */
654b5a892a1SMatthew G. Knepley   };
655ef8b56bfSJed Brown   static const PetscInt pyrArr[8*5*2] = {
656b5a892a1SMatthew G. Knepley     0, -4,  2, -3,  1, -3,  4, -3,  3, -3, /* -4: Reflect bottom face */
657b5a892a1SMatthew G. Knepley     0, -3,  3, -3,  2, -3,  1, -3,  4, -3, /* -3: Reflect bottom face */
658b5a892a1SMatthew G. Knepley     0, -2,  4, -3,  3, -3,  2, -3,  1, -3, /* -2: Reflect bottom face */
659b5a892a1SMatthew G. Knepley     0, -1,  1, -3,  4, -3,  3, -3,  2, -3, /* -1: Reflect bottom face */
660b5a892a1SMatthew G. Knepley     0,  0,  1,  0,  2,  0,  3,  0,  4,  0, /*  0: identity */
661b5a892a1SMatthew G. Knepley     0,  1,  4,  0,  1,  0,  2,  0,  3,  0, /*  1:  90 rotation about z */
662b5a892a1SMatthew G. Knepley     0,  2,  3,  0,  4,  0,  1,  0,  2,  0, /*  2: 180 rotation about z */
663b5a892a1SMatthew G. Knepley     0,  3,  2,  0,  3,  0,  4,  0,  1,  0, /*  3: 270 rotation about z */
664b5a892a1SMatthew G. Knepley   };
665b5a892a1SMatthew G. Knepley   switch (ct) {
666b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_POINT:              return pntArr;
667b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_SEGMENT:            return &segArr[(o+1)*2*2];
668b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_POINT_PRISM_TENSOR: return &segArr[(o+1)*2*2];
669b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:           return &triArr[(o+3)*3*2];
670b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:      return &quadArr[(o+4)*4*2];
671b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return &tsegArr[(o+2)*4*2];
672b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TETRAHEDRON:        return &tetArr[(o+12)*4*2];
673b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_HEXAHEDRON:         return &hexArr[(o+24)*6*2];
674b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM:          return &tripArr[(o+6)*5*2];
675b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return &ttriArr[(o+6)*5*2];
676b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return &tquadArr[(o+8)*6*2];
677b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_PYRAMID:            return &pyrArr[(o+4)*5*2];
678b5a892a1SMatthew G. Knepley     default: return NULL;
679b5a892a1SMatthew G. Knepley   }
680b5a892a1SMatthew G. Knepley }
681b5a892a1SMatthew G. Knepley 
682b5a892a1SMatthew G. Knepley /* A vertex arrangment is a vertex order */
6839fbee547SJacob Faibussowitsch static inline const PetscInt *DMPolytopeTypeGetVertexArrangment(DMPolytopeType ct, PetscInt o)
684b5a892a1SMatthew G. Knepley {
685ef8b56bfSJed Brown   static const PetscInt pntVerts[1]    = {0};
686ef8b56bfSJed Brown   static const PetscInt segVerts[2*2]  = {
687b5a892a1SMatthew G. Knepley     1, 0,
688b5a892a1SMatthew G. Knepley     0, 1};
689ef8b56bfSJed Brown   static const PetscInt triVerts[6*3]  = {
690b5a892a1SMatthew G. Knepley     1, 0, 2,
691b5a892a1SMatthew G. Knepley     0, 2, 1,
692b5a892a1SMatthew G. Knepley     2, 1, 0,
693b5a892a1SMatthew G. Knepley     0, 1, 2,
694b5a892a1SMatthew G. Knepley     1, 2, 0,
695b5a892a1SMatthew G. Knepley     2, 0, 1};
696ef8b56bfSJed Brown   static const PetscInt quadVerts[8*4]  = {
697b5a892a1SMatthew G. Knepley     2, 1, 0, 3,
698b5a892a1SMatthew G. Knepley     1, 0, 3, 2,
699b5a892a1SMatthew G. Knepley     0, 3, 2, 1,
700b5a892a1SMatthew G. Knepley     3, 2, 1, 0,
701b5a892a1SMatthew G. Knepley     0, 1, 2, 3,
702b5a892a1SMatthew G. Knepley     1, 2, 3, 0,
703b5a892a1SMatthew G. Knepley     2, 3, 0, 1,
704b5a892a1SMatthew G. Knepley     3, 0, 1, 2};
705ef8b56bfSJed Brown   static const PetscInt tsegVerts[4*4]  = {
706b5a892a1SMatthew G. Knepley     3, 2, 1, 0,
707b5a892a1SMatthew G. Knepley     1, 0, 3, 2,
708b5a892a1SMatthew G. Knepley     0, 1, 2, 3,
709b5a892a1SMatthew G. Knepley     2, 3, 0, 1};
710ef8b56bfSJed Brown   static const PetscInt tetVerts[24*4] = {
711b5a892a1SMatthew G. Knepley     2, 3, 1, 0, /* -12: (1324)   p22 */
712b5a892a1SMatthew G. Knepley     3, 1, 2, 0, /* -11: (14)     p21 */
713b5a892a1SMatthew G. Knepley     1, 2, 3, 0, /* -10: (1234)   p18 */
714b5a892a1SMatthew G. Knepley     3, 2, 0, 1, /*  -9: (1423)   p17 */
715b5a892a1SMatthew G. Knepley     2, 0, 3, 1, /*  -8: (1342)   p13 */
716b5a892a1SMatthew G. Knepley     0, 3, 2, 1, /*  -7: (24)     p14 */
717b5a892a1SMatthew G. Knepley     0, 1, 3, 2, /*  -6: (34)     p6  */
718b5a892a1SMatthew G. Knepley     1, 3, 0, 2, /*  -5: (1243)   p10 */
719b5a892a1SMatthew G. Knepley     3, 0, 1, 2, /*  -4: (1432    p9  */
720b5a892a1SMatthew G. Knepley     1, 0, 2, 3, /*  -3: (12)     p1  */
721b5a892a1SMatthew G. Knepley     0, 2, 1, 3, /*  -2: (23)     p2  */
722b5a892a1SMatthew G. Knepley     2, 1, 0, 3, /*  -1: (13)     p5  */
723b5a892a1SMatthew G. Knepley     0, 1, 2, 3, /*   0: ()       p0  */
724b5a892a1SMatthew G. Knepley     1, 2, 0, 3, /*   1: (123)    p4  */
725b5a892a1SMatthew G. Knepley     2, 0, 1, 3, /*   2: (132)    p3  */
726b5a892a1SMatthew G. Knepley     1, 0, 3, 2, /*   3: (12)(34) p7  */
727b5a892a1SMatthew G. Knepley     0, 3, 1, 2, /*   4: (243)    p8  */
728b5a892a1SMatthew G. Knepley     3, 1, 0, 2, /*   5: (143)    p11 */
729b5a892a1SMatthew G. Knepley     2, 3, 0, 1, /*   6: (13)(24) p16 */
730b5a892a1SMatthew G. Knepley     3, 0, 2, 1, /*   7: (142)    p15 */
731b5a892a1SMatthew G. Knepley     0, 2, 3, 1, /*   8: (234)    p12 */
732b5a892a1SMatthew G. Knepley     3, 2, 1, 0, /*   9: (14)(23) p23 */
733b5a892a1SMatthew G. Knepley     2, 1, 3, 0, /*  10: (134)    p19 */
734b5a892a1SMatthew G. Knepley     1, 3, 2, 0  /*  11: (124)    p20 */};
735ef8b56bfSJed Brown   static const PetscInt hexVerts[48*8] = {
736b5a892a1SMatthew G. Knepley     3,  0,  4,  5,  2,  6,  7,  1, /* -24: reflected 23 */
737b5a892a1SMatthew G. Knepley     3,  5,  6,  2,  0,  1,  7,  4, /* -23: reflected 22 */
738b5a892a1SMatthew G. Knepley     4,  0,  1,  7,  5,  6,  2,  3, /* -22: reflected 21 */
739b5a892a1SMatthew G. Knepley     6,  7,  1,  2,  5,  3,  0,  4, /* -21: reflected 20 */
740b5a892a1SMatthew G. Knepley     1,  2,  6,  7,  0,  4,  5,  3, /* -20: reflected 19 */
741b5a892a1SMatthew G. Knepley     6,  2,  3,  5,  7,  4,  0,  1, /* -19: reflected 18 */
742b5a892a1SMatthew G. Knepley     4,  5,  3,  0,  7,  1,  2,  6, /* -18: reflected 17 */
743b5a892a1SMatthew G. Knepley     1,  7,  4,  0,  2,  3,  5,  6, /* -17: reflected 16 */
744b5a892a1SMatthew G. Knepley     2,  3,  5,  6,  1,  7,  4,  0, /* -16: reflected 15 */
745b5a892a1SMatthew G. Knepley     7,  4,  0,  1,  6,  2,  3,  5, /* -15: reflected 14 */
746b5a892a1SMatthew G. Knepley     7,  1,  2,  6,  4,  5,  3,  0, /* -14: reflected 13 */
747b5a892a1SMatthew G. Knepley     0,  4,  5,  3,  1,  2,  6,  7, /* -13: reflected 12 */
748b5a892a1SMatthew G. Knepley     5,  4,  7,  6,  3,  2,  1,  0, /* -12: reflected 11 */
749b5a892a1SMatthew G. Knepley     7,  6,  5,  4,  1,  0,  3,  2, /* -11: reflected 10 */
750b5a892a1SMatthew G. Knepley     0,  1,  7,  4,  3,  5,  6,  2, /* -10: reflected  9 */
751b5a892a1SMatthew G. Knepley     4,  7,  6,  5,  0,  3,  2,  1, /*  -9: reflected  8 */
752b5a892a1SMatthew G. Knepley     5,  6,  2,  3,  4,  0,  1,  7, /*  -8: reflected  7 */
753b5a892a1SMatthew G. Knepley     2,  6,  7,  1,  3,  0,  4,  5, /*  -7: reflected  6 */
754b5a892a1SMatthew G. Knepley     6,  5,  4,  7,  2,  1,  0,  3, /*  -6: reflected  5 */
755b5a892a1SMatthew G. Knepley     5,  3,  0,  4,  6,  7,  1,  2, /*  -5: reflected  4 */
756b5a892a1SMatthew G. Knepley     2,  1,  0,  3,  6,  5,  4,  7, /*  -4: reflected  3 */
757b5a892a1SMatthew G. Knepley     1,  0,  3,  2,  7,  6,  5,  4, /*  -3: reflected  2 */
758b5a892a1SMatthew G. Knepley     0,  3,  2,  1,  4,  7,  6,  5, /*  -2: reflected  1 */
759b5a892a1SMatthew G. Knepley     3,  2,  1,  0,  5,  4,  7,  6, /*  -1: reflected  0 */
760b5a892a1SMatthew G. Knepley     0,  1,  2,  3,  4,  5,  6,  7, /*   0: identity */
761b5a892a1SMatthew G. Knepley     1,  2,  3,  0,  7,  4,  5,  6, /*   1: 90  rotation about z */
762b5a892a1SMatthew G. Knepley     2,  3,  0,  1,  6,  7,  4,  5, /*   2: 180 rotation about z */
763b5a892a1SMatthew G. Knepley     3,  0,  1,  2,  5,  6,  7,  4, /*   3: 270 rotation about z */
764b5a892a1SMatthew G. Knepley     4,  0,  3,  5,  7,  6,  2,  1, /*   4: 90  rotation about x */
765b5a892a1SMatthew G. Knepley     7,  4,  5,  6,  1,  2,  3,  0, /*   5: 180 rotation about x */
766b5a892a1SMatthew G. Knepley     1,  7,  6,  2,  0,  3,  5,  4, /*   6: 270 rotation about x */
767b5a892a1SMatthew G. Knepley     3,  2,  6,  5,  0,  4,  7,  1, /*   7: 90  rotation about y */
768b5a892a1SMatthew G. Knepley     5,  6,  7,  4,  3,  0,  1,  2, /*   8: 180 rotation about y */
769b5a892a1SMatthew G. Knepley     4,  7,  1,  0,  5,  3,  2,  6, /*   9: 270 rotation about y */
770b5a892a1SMatthew G. Knepley     4,  5,  6,  7,  0,  1,  2,  3, /*  10: 180 rotation about x+y */
771b5a892a1SMatthew G. Knepley     6,  7,  4,  5,  2,  3,  0,  1, /*  11: 180 rotation about x-y */
772b5a892a1SMatthew G. Knepley     3,  5,  4,  0,  2,  1,  7,  6, /*  12: 180 rotation about y+z */
773b5a892a1SMatthew G. Knepley     6,  2,  1,  7,  5,  4,  0,  3, /*  13: 180 rotation about y-z */
774b5a892a1SMatthew G. Knepley     1,  0,  4,  7,  2,  6,  5,  3, /*  14: 180 rotation about z+x */
775b5a892a1SMatthew G. Knepley     6,  5,  3,  2,  7,  1,  0,  4, /*  15: 180 rotation about z-x */
776b5a892a1SMatthew G. Knepley     0,  4,  7,  1,  3,  2,  6,  5, /*  16: 120 rotation about x+y+z (v0v6) */
777b5a892a1SMatthew G. Knepley     0,  3,  5,  4,  1,  7,  6,  2, /*  17: 240 rotation about x+y+z (v0v6) */
778b5a892a1SMatthew G. Knepley     5,  3,  2,  6,  4,  7,  1,  0, /*  18: 120 rotation about x+y-z (v4v2) */
779b5a892a1SMatthew G. Knepley     7,  6,  2,  1,  4,  0,  3,  5, /*  19: 240 rotation about x+y-z (v4v2) */
780b5a892a1SMatthew G. Knepley     2,  1,  7,  6,  3,  5,  4,  0, /*  20: 120 rotation about x-y+z (v1v5) */
781b5a892a1SMatthew G. Knepley     7,  1,  0,  4,  6,  5,  3,  2, /*  21: 240 rotation about x-y+z (v1v5) */
782b5a892a1SMatthew G. Knepley     2,  6,  5,  3,  1,  0,  4,  7, /*  22: 120 rotation about x-y-z (v7v3) */
783b5a892a1SMatthew G. Knepley     5,  4,  0,  3,  6,  2,  1,  7, /*  23: 240 rotation about x-y-z (v7v3) */
784b5a892a1SMatthew G. Knepley   };
785ef8b56bfSJed Brown   static const PetscInt tripVerts[12*6] = {
786b5a892a1SMatthew G. Knepley     4,  3,  5,  2,  1,  0, /* -6: reflect bottom and top */
787b5a892a1SMatthew G. Knepley     5,  4,  3,  1,  0,  2, /* -5: reflect bottom and top */
788b5a892a1SMatthew G. Knepley     3,  5,  4,  0,  2,  1, /* -4: reflect bottom and top */
789b5a892a1SMatthew G. Knepley     1,  0,  2,  5,  4,  3, /* -3: reflect bottom and top */
790b5a892a1SMatthew G. Knepley     0,  2,  1,  3,  5,  4, /* -2: reflect bottom and top */
791b5a892a1SMatthew G. Knepley     2,  1,  0,  4,  3,  5, /* -1: reflect bottom and top */
792b5a892a1SMatthew G. Knepley     0,  1,  2,  3,  4,  5, /*  0: identity */
793b5a892a1SMatthew G. Knepley     1,  2,  0,  5,  3,  4, /*  1: 120 rotation about z */
794b5a892a1SMatthew G. Knepley     2,  0,  1,  4,  5,  3, /*  2: 240 rotation about z */
795b5a892a1SMatthew G. Knepley     4,  5,  3,  2,  0,  1, /*  3: 180 rotation about y of 0 */
796b5a892a1SMatthew G. Knepley     3,  4,  5,  0,  1,  2, /*  4: 180 rotation about y of 1 */
797b5a892a1SMatthew G. Knepley     5,  3,  4,  1,  2,  0, /*  5: 180 rotation about y of 2 */
798b5a892a1SMatthew G. Knepley   };
799ef8b56bfSJed Brown   static const PetscInt ttriVerts[12*6] = {
800b5a892a1SMatthew G. Knepley     4,  3,  5,  1,  0,  2, /* -6: r b a^2 */
801b5a892a1SMatthew G. Knepley     3,  5,  4,  0,  2,  1, /* -5: r b a */
802b5a892a1SMatthew G. Knepley     5,  4,  3,  2,  1,  0, /* -4: r b */
803b5a892a1SMatthew G. Knepley     1,  0,  2,  4,  3,  5, /* -3: r a^2 */
804b5a892a1SMatthew G. Knepley     0,  2,  1,  3,  5,  4, /* -2: r a */
805b5a892a1SMatthew G. Knepley     2,  1,  0,  5,  4,  3, /* -1: r */
806b5a892a1SMatthew G. Knepley     0,  1,  2,  3,  4,  5, /*  0: identity */
807b5a892a1SMatthew G. Knepley     1,  2,  0,  4,  5,  3, /*  1: a */
808b5a892a1SMatthew G. Knepley     2,  0,  1,  5,  3,  4, /*  2: a^2 */
809b5a892a1SMatthew G. Knepley     3,  4,  5,  0,  1,  2, /*  3: b */
810b5a892a1SMatthew G. Knepley     4,  5,  3,  1,  2,  0, /*  4: b a */
811b5a892a1SMatthew G. Knepley     5,  3,  4,  2,  0,  1, /*  5: b a^2 */
812b5a892a1SMatthew G. Knepley   };
813b5a892a1SMatthew G. Knepley   /* a: rotate 90 about z
814b5a892a1SMatthew G. Knepley      b: swap top and bottom segments
815b5a892a1SMatthew G. Knepley      r: reflect */
816ef8b56bfSJed Brown   static const PetscInt tquadVerts[16*8] = {
817b5a892a1SMatthew G. Knepley     6,  5,  4,  7,  2,  1,  0,  3, /* -8: r b a^3 */
818b5a892a1SMatthew G. Knepley     5,  4,  7,  6,  1,  0,  3,  2, /* -7: r b a^2 */
819b5a892a1SMatthew G. Knepley     4,  7,  6,  5,  0,  3,  2,  1, /* -6: r b a */
820b5a892a1SMatthew G. Knepley     7,  6,  5,  4,  3,  2,  1,  0, /* -5: r b */
821b5a892a1SMatthew G. Knepley     2,  1,  0,  3,  6,  5,  4,  7, /* -4: r a^3 */
822b5a892a1SMatthew G. Knepley     1,  0,  3,  2,  5,  4,  7,  6, /* -3: r a^2 */
823b5a892a1SMatthew G. Knepley     0,  3,  2,  1,  4,  7,  6,  5, /* -2: r a */
824b5a892a1SMatthew G. Knepley     3,  2,  1,  0,  7,  6,  5,  4, /* -1: r */
825b5a892a1SMatthew G. Knepley     0,  1,  2,  3,  4,  5,  6,  7, /*  0: identity */
826b5a892a1SMatthew G. Knepley     1,  2,  3,  0,  5,  6,  7,  4, /*  1: a */
827b5a892a1SMatthew G. Knepley     2,  3,  0,  1,  6,  7,  4,  5, /*  2: a^2 */
828b5a892a1SMatthew G. Knepley     3,  0,  1,  2,  7,  4,  5,  6, /*  3: a^3 */
829b5a892a1SMatthew G. Knepley     4,  5,  6,  7,  0,  1,  2,  3, /*  4: b */
830b5a892a1SMatthew G. Knepley     5,  6,  7,  4,  1,  2,  3,  0, /*  5: b a */
831b5a892a1SMatthew G. Knepley     6,  7,  4,  5,  2,  3,  0,  1, /*  6: b a^2 */
832b5a892a1SMatthew G. Knepley     7,  4,  5,  6,  3,  0,  1,  2, /*  7: b a^3 */
833b5a892a1SMatthew G. Knepley   };
834ef8b56bfSJed Brown   static const PetscInt pyrVerts[8*5] = {
835b5a892a1SMatthew G. Knepley     2,  1,  0,  3,  4, /* -4: Reflect bottom face */
836b5a892a1SMatthew G. Knepley     1,  0,  3,  2,  4, /* -3: Reflect bottom face */
837b5a892a1SMatthew G. Knepley     0,  3,  2,  1,  4, /* -2: Reflect bottom face */
838b5a892a1SMatthew G. Knepley     3,  2,  1,  0,  4, /* -1: Reflect bottom face */
839b5a892a1SMatthew G. Knepley     0,  1,  2,  3,  4, /*  0: identity */
840b5a892a1SMatthew G. Knepley     1,  2,  3,  0,  4, /*  1:  90 rotation about z */
841b5a892a1SMatthew G. Knepley     2,  3,  0,  1,  4, /*  2: 180 rotation about z */
842b5a892a1SMatthew G. Knepley     3,  0,  1,  2,  4, /*  3: 270 rotation about z */
843b5a892a1SMatthew G. Knepley   };
844b5a892a1SMatthew G. Knepley   switch (ct) {
845b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_POINT:              return pntVerts;
846b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_SEGMENT:            return &segVerts[(o+1)*2];
847b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_POINT_PRISM_TENSOR: return &segVerts[(o+1)*2];
848b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:           return &triVerts[(o+3)*3];
849b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:      return &quadVerts[(o+4)*4];
850b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return &tsegVerts[(o+2)*4];
851b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TETRAHEDRON:        return &tetVerts[(o+12)*4];
852b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_HEXAHEDRON:         return &hexVerts[(o+24)*8];
853b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM:          return &tripVerts[(o+6)*6];
854b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return &ttriVerts[(o+6)*6];
855b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return &tquadVerts[(o+8)*8];
856b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_PYRAMID:            return &pyrVerts[(o+4)*5];
857b5a892a1SMatthew G. Knepley     default: return NULL;
858b5a892a1SMatthew G. Knepley   }
859b5a892a1SMatthew G. Knepley }
860b5a892a1SMatthew G. Knepley 
861b5a892a1SMatthew G. Knepley /* This is orientation o1 acting on orientation o2 */
8629fbee547SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeComposeOrientation(DMPolytopeType ct, PetscInt o1, PetscInt o2)
863b5a892a1SMatthew G. Knepley {
864ef8b56bfSJed Brown   static const PetscInt segMult[2*2] = {
865b5a892a1SMatthew G. Knepley      0, -1,
866b5a892a1SMatthew G. Knepley     -1,  0};
867ef8b56bfSJed Brown   static const PetscInt triMult[6*6] = {
868b5a892a1SMatthew G. Knepley      0,  2,  1, -3, -1, -2,
869b5a892a1SMatthew G. Knepley      1,  0,  2, -2, -3, -1,
870b5a892a1SMatthew G. Knepley      2,  1,  0, -1, -2, -3,
871b5a892a1SMatthew G. Knepley     -3, -2, -1,  0,  1,  2,
872b5a892a1SMatthew G. Knepley     -2, -1, -3,  1,  2,  0,
873b5a892a1SMatthew G. Knepley     -1, -3, -2,  2,  0,  1};
874ef8b56bfSJed Brown   static const PetscInt quadMult[8*8] = {
875b5a892a1SMatthew G. Knepley      0,  3,  2,  1, -4, -1, -2, -3,
876b5a892a1SMatthew G. Knepley      1,  0,  3,  2, -3, -4, -1, -2,
877b5a892a1SMatthew G. Knepley      2,  1,  0,  3, -2, -3, -4, -1,
878b5a892a1SMatthew G. Knepley      3,  2,  1,  0, -1, -2, -3, -4,
879b5a892a1SMatthew G. Knepley     -4, -3, -2, -1,  0,  1,  2,  3,
880b5a892a1SMatthew G. Knepley     -3, -2, -1, -4,  1,  2,  3,  0,
881b5a892a1SMatthew G. Knepley     -2, -1, -4, -3,  2,  3,  0,  1,
882b5a892a1SMatthew G. Knepley     -1, -4, -3, -2,  3,  0,  1,  2};
883ef8b56bfSJed Brown   static const PetscInt tsegMult[4*4] = {
884b5a892a1SMatthew G. Knepley      0,  1, -2, -1,
885b5a892a1SMatthew G. Knepley      1,  0, -1, -2,
886b5a892a1SMatthew G. Knepley     -2, -1,  0,  1,
887b5a892a1SMatthew G. Knepley     -1, -2,  1,  0};
888ef8b56bfSJed Brown   static const PetscInt tetMult[24*24] = {
889b5a892a1SMatthew G. Knepley     3, 2, 7, 0, 5, 10, 9, 8, 1, 6, 11, 4, -12, -7, -5, -9, -10, -2, -6, -1, -11, -3, -4, -8,
890b5a892a1SMatthew G. Knepley     4, 0, 8, 1, 3, 11, 10, 6, 2, 7, 9, 5, -11, -9, -4, -8, -12, -1, -5, -3, -10, -2, -6, -7,
891b5a892a1SMatthew G. Knepley     5, 1, 6, 2, 4, 9, 11, 7, 0, 8, 10, 3, -10, -8, -6, -7, -11, -3, -4, -2, -12, -1, -5, -9,
892b5a892a1SMatthew G. Knepley     0, 8, 4, 3, 11, 1, 6, 2, 10, 9, 5, 7, -9, -4, -11, -12, -1, -8, -3, -10, -5, -6, -7, -2,
893b5a892a1SMatthew G. Knepley     1, 6, 5, 4, 9, 2, 7, 0, 11, 10, 3, 8, -8, -6, -10, -11, -3, -7, -2, -12, -4, -5, -9, -1,
894b5a892a1SMatthew G. Knepley     2, 7, 3, 5, 10, 0, 8, 1, 9, 11, 4, 6, -7, -5, -12, -10, -2, -9, -1, -11, -6, -4, -8, -3,
895b5a892a1SMatthew G. Knepley     6, 5, 1, 9, 2, 4, 0, 11, 7, 3, 8, 10, -6, -10, -8, -3, -7, -11, -12, -4, -2, -9, -1, -5,
896b5a892a1SMatthew G. Knepley     7, 3, 2, 10, 0, 5, 1, 9, 8, 4, 6, 11, -5, -12, -7, -2, -9, -10, -11, -6, -1, -8, -3, -4,
897b5a892a1SMatthew G. Knepley     8, 4, 0, 11, 1, 3, 2, 10, 6, 5, 7, 9, -4, -11, -9, -1, -8, -12, -10, -5, -3, -7, -2, -6,
898b5a892a1SMatthew G. Knepley     9, 11, 10, 6, 8, 7, 3, 5, 4, 0, 2, 1, -3, -1, -2, -6, -4, -5, -9, -7, -8, -12, -10, -11,
899b5a892a1SMatthew G. Knepley     10, 9, 11, 7, 6, 8, 4, 3, 5, 1, 0, 2, -2, -3, -1, -5, -6, -4, -8, -9, -7, -11, -12, -10,
900b5a892a1SMatthew G. Knepley     11, 10, 9, 8, 7, 6, 5, 4, 3, 2, 1, 0, -1, -2, -3, -4, -5, -6, -7, -8, -9, -10, -11, -12,
901b5a892a1SMatthew G. Knepley     -12, -11, -10, -9, -8, -7, -6, -5, -4, -3, -2, -1, 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11,
902b5a892a1SMatthew G. Knepley     -11, -10, -12, -8, -7, -9, -5, -4, -6, -2, -1, -3, 1, 2, 0, 4, 5, 3, 7, 8, 6, 10, 11, 9,
903b5a892a1SMatthew G. Knepley     -10, -12, -11, -7, -9, -8, -4, -6, -5, -1, -3, -2, 2, 0, 1, 5, 3, 4, 8, 6, 7, 11, 9, 10,
904b5a892a1SMatthew G. Knepley     -9, -5, -1, -12, -2, -4, -3, -11, -7, -6, -8, -10, 3, 10, 8, 0, 7, 11, 9, 4, 2, 6, 1, 5,
905b5a892a1SMatthew G. Knepley     -8, -4, -3, -11, -1, -6, -2, -10, -9, -5, -7, -12, 4, 11, 6, 1, 8, 9, 10, 5, 0, 7, 2, 3,
906b5a892a1SMatthew G. Knepley     -7, -6, -2, -10, -3, -5, -1, -12, -8, -4, -9, -11, 5, 9, 7, 2, 6, 10, 11, 3, 1, 8, 0, 4,
907b5a892a1SMatthew G. Knepley     -3, -8, -4, -6, -11, -1, -9, -2, -10, -12, -5, -7, 6, 4, 11, 9, 1, 8, 0, 10, 5, 3, 7, 2,
908b5a892a1SMatthew G. Knepley     -2, -7, -6, -5, -10, -3, -8, -1, -12, -11, -4, -9, 7, 5, 9, 10, 2, 6, 1, 11, 3, 4, 8, 0,
909b5a892a1SMatthew G. Knepley     -1, -9, -5, -4, -12, -2, -7, -3, -11, -10, -6, -8, 8, 3, 10, 11, 0, 7, 2, 9, 4, 5, 6, 1,
910b5a892a1SMatthew G. Knepley     -6, -2, -7, -3, -5, -10, -12, -8, -1, -9, -11, -4, 9, 7, 5, 6, 10, 2, 3, 1, 11, 0, 4, 8,
911b5a892a1SMatthew G. Knepley     -5, -1, -9, -2, -4, -12, -11, -7, -3, -8, -10, -6, 10, 8, 3, 7, 11, 0, 4, 2, 9, 1, 5, 6,
912b5a892a1SMatthew G. Knepley     -4, -3, -8, -1, -6, -11, -10, -9, -2, -7, -12, -5, 11, 6, 4, 8, 9, 1, 5, 0, 10, 2, 3, 7,
913b5a892a1SMatthew G. Knepley     };
914ef8b56bfSJed Brown   static const PetscInt hexMult[48*48] = {
915b5a892a1SMatthew 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,
916b5a892a1SMatthew 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,
917b5a892a1SMatthew 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,
918b5a892a1SMatthew 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,
919b5a892a1SMatthew 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,
920b5a892a1SMatthew 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,
921b5a892a1SMatthew 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,
922b5a892a1SMatthew 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,
923b5a892a1SMatthew 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,
924b5a892a1SMatthew 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,
925b5a892a1SMatthew 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,
926b5a892a1SMatthew 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,
927b5a892a1SMatthew 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,
928b5a892a1SMatthew 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,
929b5a892a1SMatthew 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,
930b5a892a1SMatthew 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,
931b5a892a1SMatthew 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,
932b5a892a1SMatthew 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,
933b5a892a1SMatthew 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,
934b5a892a1SMatthew 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,
935b5a892a1SMatthew 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,
936b5a892a1SMatthew 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,
937b5a892a1SMatthew 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,
938b5a892a1SMatthew 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,
939b5a892a1SMatthew 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,
940b5a892a1SMatthew 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,
941b5a892a1SMatthew 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,
942b5a892a1SMatthew 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,
943b5a892a1SMatthew 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,
944b5a892a1SMatthew 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,
945b5a892a1SMatthew 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,
946b5a892a1SMatthew 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,
947b5a892a1SMatthew 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,
948b5a892a1SMatthew 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,
949b5a892a1SMatthew 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,
950b5a892a1SMatthew 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,
951b5a892a1SMatthew 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,
952b5a892a1SMatthew 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,
953b5a892a1SMatthew 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,
954b5a892a1SMatthew 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,
955b5a892a1SMatthew 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,
956b5a892a1SMatthew 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,
957b5a892a1SMatthew 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,
958b5a892a1SMatthew 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,
959b5a892a1SMatthew 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,
960b5a892a1SMatthew 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,
961b5a892a1SMatthew 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,
962b5a892a1SMatthew 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,
963b5a892a1SMatthew G. Knepley     };
964ef8b56bfSJed Brown   static const PetscInt tripMult[12*12] = {
965b5a892a1SMatthew G. Knepley     1, 0, 2, 3, 5, 4, -6, -4, -5, -2, -3, -1,
966b5a892a1SMatthew G. Knepley     0, 2, 1, 4, 3, 5, -5, -6, -4, -3, -1, -2,
967b5a892a1SMatthew G. Knepley     2, 1, 0, 5, 4, 3, -4, -5, -6, -1, -2, -3,
968b5a892a1SMatthew G. Knepley     4, 3, 5, 0, 2, 1, -3, -1, -2, -5, -6, -4,
969b5a892a1SMatthew G. Knepley     3, 5, 4, 1, 0, 2, -2, -3, -1, -6, -4, -5,
970b5a892a1SMatthew G. Knepley     5, 4, 3, 2, 1, 0, -1, -2, -3, -4, -5, -6,
971b5a892a1SMatthew G. Knepley     -6, -5, -4, -3, -2, -1, 0, 1, 2, 3, 4, 5,
972b5a892a1SMatthew G. Knepley     -4, -6, -5, -2, -1, -3, 1, 2, 0, 5, 3, 4,
973b5a892a1SMatthew G. Knepley     -5, -4, -6, -1, -3, -2, 2, 0, 1, 4, 5, 3,
974b5a892a1SMatthew G. Knepley     -3, -2, -1, -6, -5, -4, 3, 4, 5, 0, 1, 2,
975b5a892a1SMatthew G. Knepley     -1, -3, -2, -5, -4, -6, 4, 5, 3, 2, 0, 1,
976b5a892a1SMatthew G. Knepley     -2, -1, -3, -4, -6, -5, 5, 3, 4, 1, 2, 0,
977b5a892a1SMatthew G. Knepley   };
978ef8b56bfSJed Brown   static const PetscInt ttriMult[12*12] = {
979b5a892a1SMatthew G. Knepley     0, 2, 1, 3, 5, 4, -6, -4, -5, -3, -1, -2,
980b5a892a1SMatthew G. Knepley     1, 0, 2, 4, 3, 5, -5, -6, -4, -2, -3, -1,
981b5a892a1SMatthew G. Knepley     2, 1, 0, 5, 4, 3, -4, -5, -6, -1, -2, -3,
982b5a892a1SMatthew G. Knepley     3, 5, 4, 0, 2, 1, -3, -1, -2, -6, -4, -5,
983b5a892a1SMatthew G. Knepley     4, 3, 5, 1, 0, 2, -2, -3, -1, -5, -6, -4,
984b5a892a1SMatthew G. Knepley     5, 4, 3, 2, 1, 0, -1, -2, -3, -4, -5, -6,
985b5a892a1SMatthew G. Knepley     -6, -5, -4, -3, -2, -1, 0, 1, 2, 3, 4, 5,
986b5a892a1SMatthew G. Knepley     -5, -4, -6, -2, -1, -3, 1, 2, 0, 4, 5, 3,
987b5a892a1SMatthew G. Knepley     -4, -6, -5, -1, -3, -2, 2, 0, 1, 5, 3, 4,
988b5a892a1SMatthew G. Knepley     -3, -2, -1, -6, -5, -4, 3, 4, 5, 0, 1, 2,
989b5a892a1SMatthew G. Knepley     -2, -1, -3, -5, -4, -6, 4, 5, 3, 1, 2, 0,
990b5a892a1SMatthew G. Knepley     -1, -3, -2, -4, -6, -5, 5, 3, 4, 2, 0, 1,
991b5a892a1SMatthew G. Knepley   };
992ef8b56bfSJed Brown   static const PetscInt tquadMult[16*16] = {
993b5a892a1SMatthew G. Knepley     0, 3, 2, 1, 4, 7, 6, 5, -8, -5, -6, -7, -4, -1, -2, -3,
994b5a892a1SMatthew G. Knepley     1, 0, 3, 2, 5, 4, 7, 6, -7, -8, -5, -6, -3, -4, -1, -2,
995b5a892a1SMatthew G. Knepley     2, 1, 0, 3, 6, 5, 4, 7, -6, -7, -8, -5, -2, -3, -4, -1,
996b5a892a1SMatthew G. Knepley     3, 2, 1, 0, 7, 6, 5, 4, -5, -6, -7, -8, -1, -2, -3, -4,
997b5a892a1SMatthew G. Knepley     4, 7, 6, 5, 0, 3, 2, 1, -4, -1, -2, -3, -8, -5, -6, -7,
998b5a892a1SMatthew G. Knepley     5, 4, 7, 6, 1, 0, 3, 2, -3, -4, -1, -2, -7, -8, -5, -6,
999b5a892a1SMatthew G. Knepley     6, 5, 4, 7, 2, 1, 0, 3, -2, -3, -4, -1, -6, -7, -8, -5,
1000b5a892a1SMatthew G. Knepley     7, 6, 5, 4, 3, 2, 1, 0, -1, -2, -3, -4, -5, -6, -7, -8,
1001b5a892a1SMatthew G. Knepley     -8, -7, -6, -5, -4, -3, -2, -1, 0, 1, 2, 3, 4, 5, 6, 7,
1002b5a892a1SMatthew G. Knepley     -7, -6, -5, -8, -3, -2, -1, -4, 1, 2, 3, 0, 5, 6, 7, 4,
1003b5a892a1SMatthew G. Knepley     -6, -5, -8, -7, -2, -1, -4, -3, 2, 3, 0, 1, 6, 7, 4, 5,
1004b5a892a1SMatthew G. Knepley     -5, -8, -7, -6, -1, -4, -3, -2, 3, 0, 1, 2, 7, 4, 5, 6,
1005b5a892a1SMatthew G. Knepley     -4, -3, -2, -1, -8, -7, -6, -5, 4, 5, 6, 7, 0, 1, 2, 3,
1006b5a892a1SMatthew G. Knepley     -3, -2, -1, -4, -7, -6, -5, -8, 5, 6, 7, 4, 1, 2, 3, 0,
1007b5a892a1SMatthew G. Knepley     -2, -1, -4, -3, -6, -5, -8, -7, 6, 7, 4, 5, 2, 3, 0, 1,
1008b5a892a1SMatthew G. Knepley     -1, -4, -3, -2, -5, -8, -7, -6, 7, 4, 5, 6, 3, 0, 1, 2,
1009b5a892a1SMatthew G. Knepley   };
1010ef8b56bfSJed Brown   static const PetscInt pyrMult[8*8] = {
1011b5a892a1SMatthew G. Knepley     0, 3, 2, 1, -4, -1, -2, -3,
1012b5a892a1SMatthew G. Knepley     1, 0, 3, 2, -3, -4, -1, -2,
1013b5a892a1SMatthew G. Knepley     2, 1, 0, 3, -2, -3, -4, -1,
1014b5a892a1SMatthew G. Knepley     3, 2, 1, 0, -1, -2, -3, -4,
1015b5a892a1SMatthew G. Knepley     -4, -3, -2, -1, 0, 1, 2, 3,
1016b5a892a1SMatthew G. Knepley     -3, -2, -1, -4, 1, 2, 3, 0,
1017b5a892a1SMatthew G. Knepley     -2, -1, -4, -3, 2, 3, 0, 1,
1018b5a892a1SMatthew G. Knepley     -1, -4, -3, -2, 3, 0, 1, 2,
1019b5a892a1SMatthew G. Knepley   };
1020b5a892a1SMatthew G. Knepley   switch (ct) {
1021b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_POINT:              return 0;
1022b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_SEGMENT:
1023b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_POINT_PRISM_TENSOR: return segMult[(o1+1)*2+o2+1];
1024b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:           return triMult[(o1+3)*6+o2+3];
1025b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:      return quadMult[(o1+4)*8+o2+4];
1026b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return tsegMult[(o1+2)*4+o2+2];
1027b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TETRAHEDRON:        return tetMult[(o1+12)*24+o2+12];
1028b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_HEXAHEDRON:         return hexMult[(o1+24)*48+o2+24];
1029b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM:          return tripMult[(o1+6)*12+o2+6];
1030b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return ttriMult[(o1+6)*12+o2+6];
1031b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return tquadMult[(o1+8)*16+o2+8];
1032b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_PYRAMID:            return pyrMult[(o1+4)*8+o2+4];
1033b5a892a1SMatthew G. Knepley     default: return 0;
1034b5a892a1SMatthew G. Knepley   }
1035b5a892a1SMatthew G. Knepley }
1036b5a892a1SMatthew G. Knepley 
1037b5a892a1SMatthew G. Knepley /* This is orientation o1 acting on orientation o2^{-1} */
10389fbee547SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeComposeOrientationInv(DMPolytopeType ct, PetscInt o1, PetscInt o2)
1039b5a892a1SMatthew G. Knepley {
1040ef8b56bfSJed Brown   static const PetscInt triInv[6]    = {-3, -2, -1, 0, 2, 1};
1041ef8b56bfSJed Brown   static const PetscInt quadInv[8]   = {-4, -3, -2, -1, 0, 3, 2, 1};
1042ef8b56bfSJed 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};
1043ef8b56bfSJed Brown   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,
1044b5a892a1SMatthew G. Knepley                                           0,   3,   2,   1,   6,   5,   4,   9,   8,   7,  10,  11,  12,  13,  14, 15, 17, 16, 19, 18, 21, 20, 23, 22};
1045ef8b56bfSJed Brown   static const PetscInt tripInv[12]  = {-5, -6, -4, -3, -2, -1, 0, 2, 1, 3, 4, 5};
1046ef8b56bfSJed Brown   static const PetscInt ttriInv[12]  = {-6, -5, -4, -3, -2, -1, 0, 2, 1, 3, 5, 4};
1047ef8b56bfSJed Brown   static const PetscInt tquadInv[16] = {-8, -7, -6, -5, -4, -3, -2, -1, 0, 3, 2, 1, 4, 7, 6, 5};
1048ef8b56bfSJed Brown   static const PetscInt pyrInv[8]    = {-4, -3, -2, -1, 0, 3, 2, 1};
1049b5a892a1SMatthew G. Knepley   switch (ct) {
1050b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_POINT:              return 0;
1051b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_SEGMENT:
1052b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_POINT_PRISM_TENSOR: return DMPolytopeTypeComposeOrientation(ct, o1, o2);
1053b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:           return DMPolytopeTypeComposeOrientation(ct, o1, triInv[o2+3]);
1054b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:      return DMPolytopeTypeComposeOrientation(ct, o1, quadInv[o2+4]);
1055b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return DMPolytopeTypeComposeOrientation(ct, o1, o2);
1056b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TETRAHEDRON:        return DMPolytopeTypeComposeOrientation(ct, o1, tetInv[o2+12]);
1057b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_HEXAHEDRON:         return DMPolytopeTypeComposeOrientation(ct, o1, hexInv[o2+24]);
1058b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM:          return DMPolytopeTypeComposeOrientation(ct, o1, tripInv[o2+6]);
1059b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return DMPolytopeTypeComposeOrientation(ct, o1, ttriInv[o2+6]);
1060b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return DMPolytopeTypeComposeOrientation(ct, o1, tquadInv[o2+8]);
1061b5a892a1SMatthew G. Knepley     case DM_POLYTOPE_PYRAMID:            return DMPolytopeTypeComposeOrientation(ct, o1, pyrInv[o2+4]);
1062b5a892a1SMatthew G. Knepley     default: return 0;
1063b5a892a1SMatthew G. Knepley   }
1064b5a892a1SMatthew G. Knepley }
1065b5a892a1SMatthew G. Knepley 
1066b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeMatchOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *);
1067b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeMatchVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *);
1068b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeGetOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *);
1069b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeGetVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *);
1070012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeInCellTest(DMPolytopeType, const PetscReal[], PetscBool *);
1071b5a892a1SMatthew G. Knepley 
1072e1589f56SBarry Smith #endif
1073