xref: /petsc/include/petscdm.h (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1e1589f56SBarry Smith /*
2e1589f56SBarry Smith       Objects to manage the interactions between the mesh data structures and the algebraic objects
3e1589f56SBarry Smith */
4a4963045SJacob Faibussowitsch #pragma once
5bb4b53efSMatthew G. Knepley #include "petscsystypes.h"
62c8e378dSBarry Smith #include <petscmat.h>
71e25c274SJed Brown #include <petscdmtypes.h>
8ce78bad3SBarry Smith #include <petscdmlabel.h>
9decb47aaSMatthew G. Knepley #include <petscfetypes.h>
102764a2aaSMatthew G. Knepley #include <petscdstypes.h>
1107218a29SMatthew G. Knepley #include <petscdt.h>
12e1589f56SBarry Smith 
13ac09b921SBarry Smith /* SUBMANSEC = DM */
14ac09b921SBarry Smith 
15607a6623SBarry Smith PETSC_EXTERN PetscErrorCode DMInitializePackage(void);
16e1589f56SBarry Smith 
17014dd563SJed Brown PETSC_EXTERN PetscClassId DM_CLASSID;
18e1589f56SBarry Smith 
19528c63e0SDave May #define DMLOCATEPOINT_POINT_NOT_FOUND -367
20528c63e0SDave May 
2176bdecfbSBarry Smith /*J
220b4b7b1cSBarry Smith    DMType - String with the name of a PETSc `DM`. These are all the `DM` provided by PETSc.
23e1589f56SBarry Smith 
24e1589f56SBarry Smith    Level: beginner
25e1589f56SBarry Smith 
260b4b7b1cSBarry Smith    Note:
270b4b7b1cSBarry Smith    These can be used with `DMSetType()` or the options database key `-dm_type` to set the specific data structures and algorithms to use with a specific `DM`.
280b4b7b1cSBarry Smith    But more commonly one calls directly a constructor for a particular `DMType` such as `DMDACreate()`
290b4b7b1cSBarry Smith 
300b4b7b1cSBarry Smith .seealso: [](ch_dmbase), `DMSetType()`, `DMCreate()`, `DM`, `DMDACreate()`
3176bdecfbSBarry Smith J*/
3219fd82e9SBarry Smith typedef const char *DMType;
33e1589f56SBarry Smith #define DMDA        "da"
34e1589f56SBarry Smith #define DMCOMPOSITE "composite"
35e1589f56SBarry Smith #define DMSLICED    "sliced"
36fe1899a2SJed Brown #define DMSHELL     "shell"
37ab7f58a0SBarry Smith #define DMPLEX      "plex"
388ac4e037SJed Brown #define DMREDUNDANT "redundant"
393a19ef87SMatthew G Knepley #define DMPATCH     "patch"
401d72bce8STim Tautges #define DMMOAB      "moab"
415f2c45f1SShri Abhyankar #define DMNETWORK   "network"
42ef51cf95SToby Isaac #define DMFOREST    "forest"
43db4d5e8cSToby Isaac #define DMP4EST     "p4est"
44db4d5e8cSToby Isaac #define DMP8EST     "p8est"
452fd35b1fSDave May #define DMSWARM     "swarm"
46d852a638SPatrick Sanan #define DMPRODUCT   "product"
47a3101111SPatrick Sanan #define DMSTAG      "stag"
48e1589f56SBarry Smith 
49bff4a2f0SMatthew G. Knepley PETSC_EXTERN const char *const       DMBoundaryTypes[];
5040967b3bSMatthew G. Knepley PETSC_EXTERN const char *const       DMBoundaryConditionTypes[];
51863027abSJed Brown PETSC_EXTERN const char *const       DMBlockingTypes[];
52140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList       DMList;
53c0517cd5SMatthew G. Knepley PETSC_EXTERN DMGeneratorFunctionList DMGenerateList;
547625649eSMatthew G. Knepley PETSC_EXTERN PetscFunctionList       DMGeomModelList;
55014dd563SJed Brown PETSC_EXTERN PetscErrorCode          DMCreate(MPI_Comm, DM *);
5638221697SMatthew G. Knepley PETSC_EXTERN PetscErrorCode          DMClone(DM, DM *);
5719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode          DMSetType(DM, DMType);
5819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode          DMGetType(DM, DMType *);
59bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode          DMRegister(const char[], PetscErrorCode (*)(DM));
60014dd563SJed Brown PETSC_EXTERN PetscErrorCode          DMRegisterDestroy(void);
61e1589f56SBarry Smith 
62014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMView(DM, PetscViewer);
63014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLoad(DM, PetscViewer);
64014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDestroy(DM *);
65014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM, Vec *);
66014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM, Vec *);
67014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM, Vec *);
68014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM, Vec *);
69014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM, Vec *);
70014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM, Vec *);
71014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
7250eeb1caSToby Isaac PETSC_EXTERN PetscErrorCode DMClearLocalVectors(DM);
73974ca4ecSStefano Zampini PETSC_EXTERN PetscErrorCode DMClearNamedGlobalVectors(DM);
74974ca4ecSStefano Zampini PETSC_EXTERN PetscErrorCode DMClearNamedLocalVectors(DM);
75e77ac854SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasNamedGlobalVector(DM, const char *, PetscBool *);
76014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM, const char *, Vec *);
77014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM, const char *, Vec *);
78e77ac854SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasNamedLocalVector(DM, const char *, PetscBool *);
792348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM, const char *, Vec *);
802348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM, const char *, Vec *);
81014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM, ISLocalToGlobalMapping *);
82014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM, PetscInt *, char ***, IS **);
83014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM, PetscInt *);
84b412c318SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateColoring(DM, ISColoringType, ISColoring *);
85b412c318SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM, Mat *);
86aa0f6e3cSJed Brown PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateSkip(DM, PetscBool);
87014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM, PetscBool);
88b06ff27eSHong Zhang PETSC_EXTERN PetscErrorCode DMSetMatrixStructureOnly(DM, PetscBool);
89863027abSJed Brown PETSC_EXTERN PetscErrorCode DMSetBlockingType(DM, DMBlockingType);
90863027abSJed Brown PETSC_EXTERN PetscErrorCode DMGetBlockingType(DM, DMBlockingType *);
91014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM, DM, Mat *, Vec *);
923ad4599aSBarry Smith PETSC_EXTERN PetscErrorCode DMCreateRestriction(DM, DM, Mat *);
936dbf9973SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMCreateInjection(DM, DM, Mat *);
94bd041c0cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateMassMatrix(DM, DM, Mat *);
958e9849d2SStefano Zampini PETSC_EXTERN PetscErrorCode DMCreateMassMatrixLumped(DM, Vec *, Vec *);
961898fd5cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateGradientMatrix(DM, DM, Mat *);
9769291d52SBarry Smith PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM, PetscInt, MPI_Datatype, void *);
9869291d52SBarry Smith PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM, PetscInt, MPI_Datatype, void *);
99014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefine(DM, MPI_Comm, DM *);
100014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsen(DM, MPI_Comm, DM *);
101a8fb8f29SToby Isaac PETSC_EXTERN PetscErrorCode DMGetCoarseDM(DM, DM *);
102a8fb8f29SToby Isaac PETSC_EXTERN PetscErrorCode DMSetCoarseDM(DM, DM);
10388bdff64SToby Isaac PETSC_EXTERN PetscErrorCode DMGetFineDM(DM, DM *);
10488bdff64SToby Isaac PETSC_EXTERN PetscErrorCode DMSetFineDM(DM, DM);
105014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM, PetscInt, DM[]);
106014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM, PetscInt, DM[]);
107*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM, PetscErrorCode (*)(DM, DM, PetscCtx), PetscErrorCode (*)(DM, Mat, Vec, Mat, DM, PetscCtx), PetscCtx);
108*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMCoarsenHookRemove(DM, PetscErrorCode (*)(DM, DM, PetscCtx), PetscErrorCode (*)(DM, Mat, Vec, Mat, DM, PetscCtx), PetscCtx);
109*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM, PetscErrorCode (*)(DM, DM, PetscCtx), PetscErrorCode (*)(DM, Mat, DM, PetscCtx), PetscCtx);
110*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMRefineHookRemove(DM, PetscErrorCode (*)(DM, DM, PetscCtx), PetscErrorCode (*)(DM, Mat, DM, PetscCtx), PetscCtx);
111014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestrict(DM, Mat, Vec, Mat, DM);
112014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMInterpolate(DM, Mat, DM);
1131f3379b2SToby Isaac PETSC_EXTERN PetscErrorCode DMInterpolateSolution(DM, DM, Mat, Vec, Vec);
114d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMExtrude(DM, PetscInt, DM *);
115014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
116fe2efc57SMark PETSC_EXTERN PetscErrorCode DMViewFromOptions(DM, PetscObject, const char[]);
117ca266f36SBarry Smith 
118c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerate(DM, const char[], PetscBool, DM *);
1199fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMGenerateRegister(const char[], PetscErrorCode (*)(DM, PetscBool, DM *), PetscErrorCode (*)(DM, PetscReal *, DM *), PetscErrorCode (*)(DM, Vec, DMLabel, DMLabel, DM *), PetscInt);
120c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerateRegisterAll(void);
121c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerateRegisterDestroy(void);
1227625649eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGeomModelRegister(const char[], PetscErrorCode (*)(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]));
1237625649eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGeomModelRegisterAll(void);
1247625649eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGeomModelRegisterDestroy(void);
125a1b0c543SToby Isaac PETSC_EXTERN PetscErrorCode DMAdaptLabel(DM, DMLabel, DM *);
1269fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMAdaptMetric(DM, Vec, DMLabel, DMLabel, DM *);
127df0b854cSToby Isaac 
128014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetUp(DM);
129014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM, DM, Mat, Vec *);
130edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMDACreateAggregates()", ) PetscErrorCode DMCreateAggregates(DM, DM, Mat *);
131*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM, PetscErrorCode (*)(DM, Vec, InsertMode, Vec, PetscCtx), PetscErrorCode (*)(DM, Vec, InsertMode, Vec, PetscCtx), PetscCtx);
132*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMLocalToGlobalHookAdd(DM, PetscErrorCode (*)(DM, Vec, InsertMode, Vec, PetscCtx), PetscErrorCode (*)(DM, Vec, InsertMode, Vec, PetscCtx), PetscCtx);
13301729b5cSPatrick Sanan PETSC_EXTERN PetscErrorCode DMGlobalToLocal(DM, Vec, InsertMode, Vec);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM, Vec, InsertMode, Vec);
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM, Vec, InsertMode, Vec);
13601729b5cSPatrick Sanan PETSC_EXTERN PetscErrorCode DMLocalToGlobal(DM, Vec, InsertMode, Vec);
137014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM, Vec, InsertMode, Vec);
138014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM, Vec, InsertMode, Vec);
139d78e899eSRichard Tran Mills PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM, Vec, InsertMode, Vec);
140d78e899eSRichard Tran Mills PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM, Vec, InsertMode, Vec);
14119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMConvert(DM, DMType, DM *);
142e1589f56SBarry Smith 
143c73cfb54SMatthew G. Knepley /* Topology support */
144c73cfb54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDimension(DM, PetscInt *);
145c73cfb54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetDimension(DM, PetscInt);
146793f3fe5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDimPoints(DM, PetscInt, PetscInt *, PetscInt *);
1478e4ac7eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetUseNatural(DM, PetscBool *);
1488e4ac7eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetUseNatural(DM, PetscBool);
149ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode DMGetNeighbors(DM, PetscInt *, const PetscMPIInt *[]);
15046e270d4SMatthew G. Knepley 
15146e270d4SMatthew G. Knepley /* Coordinate support */
1526636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM, DM *);
1531cfe2091SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM, DM);
1546858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinateDM(DM, DM *);
1556858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinateDM(DM, DM);
15646e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateDim(DM, PetscInt *);
15746e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateDim(DM, PetscInt);
158e8abe2deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM, PetscSection *);
15946e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM, PetscInt, PetscSection);
1606858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinateSection(DM, PetscSection *);
1616858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinateSection(DM, PetscInt, PetscSection);
1626636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM, Vec *);
1636636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM, Vec);
1646858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinates(DM, Vec *);
1656858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinates(DM, Vec);
16681e9a530SVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalSetUp(DM);
1676858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM, Vec *);
16881e9a530SVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalNoncollective(DM, Vec *);
1692db98f8dSVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalTuple(DM, IS, PetscSection *, Vec *);
1706636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM, Vec);
1716858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM);
1726858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocal(DM, Vec *);
1736858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM, Vec *);
1746858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinatesLocal(DM, Vec);
1756858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateField(DM, DMField *);
1766858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateField(DM, DMField);
1774c712d99Sksagiyam PETSC_EXTERN PetscErrorCode DMSetCellCoordinateField(DM, DMField);
1786858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetLocalBoundingBox(DM, PetscReal[], PetscReal[]);
1796858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetBoundingBox(DM, PetscReal[], PetscReal[]);
1804c712d99Sksagiyam PETSC_EXTERN PetscErrorCode DMSetCoordinateDisc(DM, PetscFE, PetscBool, PetscBool);
18162a38674SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocatePoints(DM, Vec, DMPointLocationType, PetscSF *);
1827625649eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSnapToGeomModel(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
1837625649eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetSnapToGeomModel(DM, const char[]);
1846858538eSMatthew G. Knepley 
1856858538eSMatthew G. Knepley /* Periodicity support */
1864fb89dddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM, const PetscReal *[], const PetscReal *[], const PetscReal *[]);
1874fb89dddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM, const PetscReal[], const PetscReal[], const PetscReal[]);
188e907e85cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate(DM, const PetscScalar[], PetscBool, PetscScalar[]);
1892e17dfb7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocalizeCoordinates(DM);
19036447a5eSToby Isaac PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalized(DM, PetscBool *);
1918f700142SStefano Zampini PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalizedLocal(DM, PetscBool *);
192a4b8866cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetSparseLocalize(DM, PetscBool *);
193a4b8866cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetSparseLocalize(DM, PetscBool);
1946636e97aSMatthew G Knepley 
1955dbd56e3SPeter Brune /* block hook interface */
196*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM, PetscErrorCode (*)(DM, DM, PetscCtx), PetscErrorCode (*)(DM, VecScatter, VecScatter, DM, PetscCtx), PetscCtx);
197*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSubDomainHookRemove(DM, PetscErrorCode (*)(DM, DM, PetscCtx), PetscErrorCode (*)(DM, VecScatter, VecScatter, DM, PetscCtx), PetscCtx);
198be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM, VecScatter, VecScatter, DM);
1995dbd56e3SPeter Brune 
200014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM, const char[]);
20131697293SDave May PETSC_EXTERN PetscErrorCode DMAppendOptionsPrefix(DM, const char[]);
20231697293SDave May PETSC_EXTERN PetscErrorCode DMGetOptionsPrefix(DM, const char *[]);
20319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetVecType(DM, VecType);
204c0dedaeaSBarry Smith PETSC_EXTERN PetscErrorCode DMGetVecType(DM, VecType *);
20519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetMatType(DM, MatType);
206c0dedaeaSBarry Smith PETSC_EXTERN PetscErrorCode DMGetMatType(DM, MatType *);
2078f1509bcSBarry Smith PETSC_EXTERN PetscErrorCode DMSetISColoringType(DM, ISColoringType);
2088f1509bcSBarry Smith PETSC_EXTERN PetscErrorCode DMGetISColoringType(DM, ISColoringType *);
209*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM, PetscCtx);
21049abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM, PetscCtxDestroyFn *);
211*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM, PetscCtxRt);
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM, PetscErrorCode (*)(DM, Vec, Vec));
213014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM, PetscBool *);
214b0ae01b7SPeter Brune PETSC_EXTERN PetscErrorCode DMHasColoring(DM, PetscBool *);
2153ad4599aSBarry Smith PETSC_EXTERN PetscErrorCode DMHasCreateRestriction(DM, PetscBool *);
216a7058e45SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMHasCreateInjection(DM, PetscBool *);
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM, Vec, Vec);
21893d92d96SBarry Smith 
21937bc7515SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, const PetscInt[], IS *, DM *);
2202adcc780SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSuperDM(DM[], PetscInt, IS **, DM *);
221dd072f5fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSectionSubDM(DM, PetscInt, const PetscInt[], const PetscInt[], const PetscInt[], IS *, DM *);
222ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateSectionSuperDM(DM[], PetscInt, IS *[], DM *);
22316621825SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM, PetscInt *, char ***, IS **, DM **);
2248d4ac253SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM, PetscInt *, char ***, IS **, IS **, DM **);
225ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM, PetscInt, DM *, VecScatter *[], VecScatter *[], VecScatter *[]);
226e7c4fc90SDmitry Karpeev 
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM, PetscInt *);
228fef3a512SBarry Smith PETSC_EXTERN PetscErrorCode DMSetRefineLevel(DM, PetscInt);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM, PetscInt *);
2309a64c4a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoarsenLevel(DM, PetscInt);
231014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);
232e1589f56SBarry Smith 
2335f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM *);
2345f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
235c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM *);
236c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
237531c7667SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringUseDM(Mat, MatFDColoring);
2385f1ad066SMatthew G Knepley 
239e1589f56SBarry Smith typedef struct NLF_DAAD *NLF;
240e1589f56SBarry Smith 
241bc2bf880SBarry Smith #define DM_FILE_CLASSID 1211221
2427da65231SMatthew G Knepley 
2437da65231SMatthew G Knepley /* FEM support */
2442ecf6ec3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPrintCellIndices(PetscInt, const char[], PetscInt, const PetscInt[]);
245014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char[], PetscInt, const PetscScalar[]);
2465962854dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPrintCellVectorReal(PetscInt, const char[], PetscInt, const PetscReal[]);
247014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char[], PetscInt, PetscInt, const PetscScalar[]);
2486113b454SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char[], PetscReal, Vec);
2497da65231SMatthew G Knepley 
2508cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *));
2518cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *));
2528cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *));
2538cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *));
254f9d4088aSMatthew G. Knepley 
255061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalSection(DM, PetscSection *);
256061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMSetLocalSection(DM, PetscSection);
257e87a4003SBarry Smith PETSC_EXTERN PetscErrorCode DMGetGlobalSection(DM, PetscSection *);
258e87a4003SBarry Smith PETSC_EXTERN PetscErrorCode DMSetGlobalSection(DM, PetscSection);
259adc21957SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSectionPermutation(DM, IS *, PetscBT *);
260adc21957SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMReorderSectionGetDefault(DM, DMReorderDefaultFlag *);
261adc21957SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMReorderSectionSetDefault(DM, DMReorderDefaultFlag);
262adc21957SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMReorderSectionGetType(DM, MatOrderingType *);
263adc21957SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMReorderSectionSetType(DM, MatOrderingType);
264d2b2dc1eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMUseTensorOrder(DM, PetscBool);
DMGetSection(DM dm,PetscSection * s)265b98a7184SJames Wright static inline PETSC_DEPRECATED_FUNCTION(3, 23, 0, "DMGetLocalSection()", ) PetscErrorCode DMGetSection(DM dm, PetscSection *s)
266b98a7184SJames Wright {
267b98a7184SJames Wright   return DMGetLocalSection(dm, s);
268b98a7184SJames Wright }
DMSetSection(DM dm,PetscSection s)269b98a7184SJames Wright static inline PETSC_DEPRECATED_FUNCTION(3, 23, 0, "DMSetLocalSection()", ) PetscErrorCode DMSetSection(DM dm, PetscSection s)
270b98a7184SJames Wright {
271b98a7184SJames Wright   return DMSetLocalSection(dm, s);
272b98a7184SJames Wright }
DMGetDefaultSection(DM dm,PetscSection * s)273edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMGetSection()", ) PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *s)
274d71ae5a4SJacob Faibussowitsch {
275b98a7184SJames Wright   return DMGetLocalSection(dm, s);
2769371c9d4SSatish Balay }
DMSetDefaultSection(DM dm,PetscSection s)277edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMSetSection()", ) PetscErrorCode DMSetDefaultSection(DM dm, PetscSection s)
278d71ae5a4SJacob Faibussowitsch {
279b98a7184SJames Wright   return DMSetLocalSection(dm, s);
2809371c9d4SSatish Balay }
DMGetDefaultGlobalSection(DM dm,PetscSection * s)281edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMGetGlobalSection()", ) PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *s)
282d71ae5a4SJacob Faibussowitsch {
2839371c9d4SSatish Balay   return DMGetGlobalSection(dm, s);
2849371c9d4SSatish Balay }
DMSetDefaultGlobalSection(DM dm,PetscSection s)285edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMSetGlobalSection()", ) PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection s)
286d71ae5a4SJacob Faibussowitsch {
2879371c9d4SSatish Balay   return DMSetGlobalSection(dm, s);
2889371c9d4SSatish Balay }
289e87a4003SBarry Smith 
2901bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMGetSectionSF(DM, PetscSF *);
2911bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMSetSectionSF(DM, PetscSF);
2921bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateSectionSF(DM, PetscSection, PetscSection);
DMGetDefaultSF(DM dm,PetscSF * s)293edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMGetSectionSF()", ) PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *s)
294d71ae5a4SJacob Faibussowitsch {
2959371c9d4SSatish Balay   return DMGetSectionSF(dm, s);
2969371c9d4SSatish Balay }
DMSetDefaultSF(DM dm,PetscSF s)297edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMSetSectionSF()", ) PetscErrorCode DMSetDefaultSF(DM dm, PetscSF s)
298d71ae5a4SJacob Faibussowitsch {
2999371c9d4SSatish Balay   return DMSetSectionSF(dm, s);
3009371c9d4SSatish Balay }
DMCreateDefaultSF(DM dm,PetscSection l,PetscSection g)301edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMCreateSectionSF()", ) PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection l, PetscSection g)
302d71ae5a4SJacob Faibussowitsch {
3039371c9d4SSatish Balay   return DMCreateSectionSF(dm, l, g);
3049371c9d4SSatish Balay }
3051bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
3061bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);
3074f37162bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNaturalSF(DM, PetscSF *);
3084f37162bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNaturalSF(DM, PetscSF);
3091bb6d2a8SBarry Smith 
31079769bd5SJed Brown PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *, Vec *);
31179769bd5SJed Brown PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat, Vec);
312e87a4003SBarry Smith 
31314f150ffSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *);
314cdb7a50dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *);
315cdb7a50dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal);
316b2033f5dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char[], PetscInt, PetscReal *);
317b2033f5dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetOutputSequenceLength(DM, PetscViewer, const char[], PetscInt *);
31814f150ffSMatthew G. Knepley 
319af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
320af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
32144a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, DMLabel *, PetscObject *);
32244a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, DMLabel, PetscObject);
32344a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAddField(DM, DMLabel, PetscObject);
324e0b68406SMatthew Knepley PETSC_EXTERN PetscErrorCode DMSetFieldAvoidTensor(DM, PetscInt, PetscBool);
325e0b68406SMatthew Knepley PETSC_EXTERN PetscErrorCode DMGetFieldAvoidTensor(DM, PetscInt, PetscBool *);
32644a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMClearFields(DM);
327bb4b53efSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyFields(DM, PetscInt, PetscInt, DM);
32834aa8a36SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAdjacency(DM, PetscInt, PetscBool *, PetscBool *);
32934aa8a36SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetAdjacency(DM, PetscInt, PetscBool, PetscBool);
330b0441da4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetBasicAdjacency(DM, PetscBool *, PetscBool *);
331b0441da4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetBasicAdjacency(DM, PetscBool, PetscBool);
332e5e52638SMatthew G. Knepley 
333e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNumDS(DM, PetscInt *);
334e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *);
33507218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellDS(DM, PetscInt, PetscDS *, PetscDS *);
33607218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetRegionDS(DM, DMLabel, IS *, PetscDS *, PetscDS *);
33707218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetRegionDS(DM, DMLabel, IS, PetscDS, PetscDS);
33807218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetRegionNumDS(DM, PetscInt, DMLabel *, IS *, PetscDS *, PetscDS *);
33907218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetRegionNumDS(DM, PetscInt, DMLabel, IS, PetscDS, PetscDS);
3401d3af9e0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMFindRegionNum(DM, PetscDS, PetscInt *);
3412df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateFEDefault(DM, PetscInt, const char[], PetscInt, PetscFE *);
342e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateDS(DM);
343e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMClearDS(DM);
344bb4b53efSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyDS(DM, PetscInt, PetscInt, DM);
345e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyDisc(DM, DM);
346f2cacb80SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMComputeExactSolution(DM, PetscReal, Vec, Vec);
3479a2a23afSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNumAuxiliaryVec(DM, PetscInt *);
348ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAuxiliaryVec(DM, DMLabel, PetscInt, PetscInt, Vec *);
349ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetAuxiliaryVec(DM, DMLabel, PetscInt, PetscInt, Vec);
350ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAuxiliaryLabels(DM, DMLabel[], PetscInt[], PetscInt[]);
3519a2a23afSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyAuxiliaryVec(DM, DM);
352e4d5475eSStefano Zampini PETSC_EXTERN PetscErrorCode DMClearAuxiliaryVec(DM);
353af122d2aSMatthew G Knepley 
3544267b1a3SMatthew G. Knepley /*MC
355ce78bad3SBarry Smith   DMInterpolationInfo - Pointer to a structure for holding information about interpolation on a mesh
3564267b1a3SMatthew G. Knepley 
3574267b1a3SMatthew G. Knepley   Synopsis:
3584267b1a3SMatthew G. Knepley     comm     - The communicator
3594267b1a3SMatthew G. Knepley     dim      - The spatial dimension of points
3604267b1a3SMatthew G. Knepley     nInput   - The number of input points
361ce78bad3SBarry Smith     points[] - The input point coordinates
362ce78bad3SBarry Smith     cells[]  - The cell containing each point
3634267b1a3SMatthew G. Knepley     n        - The number of local points
3644267b1a3SMatthew G. Knepley     coords   - The point coordinates
3654267b1a3SMatthew G. Knepley     dof      - The number of components to interpolate
3664267b1a3SMatthew G. Knepley 
36716a05f60SBarry Smith   Level: intermediate
36816a05f60SBarry Smith 
369af27ebaaSBarry Smith .seealso: [](ch_dmbase), `DM`, `DMInterpolationCreate()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`
3704267b1a3SMatthew G. Knepley M*/
371ce78bad3SBarry Smith struct _n_DMInterpolationInfo {
372e87bb0d3SMatthew G Knepley   MPI_Comm   comm;
373ce78bad3SBarry Smith   PetscInt   dim;    /* The spatial dimension of points */
374e87bb0d3SMatthew G Knepley   PetscInt   nInput; /* The number of input points */
375e87bb0d3SMatthew G Knepley   PetscReal *points; /* The input point coordinates */
376e87bb0d3SMatthew G Knepley   PetscInt  *cells;  /* The cell containing each point */
377e87bb0d3SMatthew G Knepley   PetscInt   n;      /* The number of local points */
378e87bb0d3SMatthew G Knepley   Vec        coords; /* The point coordinates */
379e87bb0d3SMatthew G Knepley   PetscInt   dof;    /* The number of components to interpolate */
380e87bb0d3SMatthew G Knepley };
381ce78bad3SBarry Smith typedef struct _n_DMInterpolationInfo *DMInterpolationInfo;
382e87bb0d3SMatthew G Knepley 
38394b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
38494b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
38594b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
38694b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
38794b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
38894b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
38952aa1562SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool, PetscBool);
39094b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
39194b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
39294b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
39394b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
39494b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);
395c58f1c22SToby Isaac 
396c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMCreateLabel(DM, const char[]);
3974bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreateLabelAtIndex(DM, PetscInt, const char[]);
398c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelValue(DM, const char[], PetscInt, PetscInt *);
399c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMSetLabelValue(DM, const char[], PetscInt, PetscInt);
400c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMClearLabelValue(DM, const char[], PetscInt, PetscInt);
401c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelSize(DM, const char[], PetscInt *);
402c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelIdIS(DM, const char[], IS *);
403c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetStratumSize(DM, const char[], PetscInt, PetscInt *);
404c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetStratumIS(DM, const char[], PetscInt, IS *);
4054de306b1SToby Isaac PETSC_EXTERN PetscErrorCode DMSetStratumIS(DM, const char[], PetscInt, IS);
406c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMClearLabelStratum(DM, const char[], PetscInt);
407c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelOutput(DM, const char[], PetscBool *);
408c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMSetLabelOutput(DM, const char[], PetscBool);
4095f15299fSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMGetFirstLabeledPoint(DM, DM, DMLabel, PetscInt, const PetscInt *, PetscInt, PetscInt *, PetscDS *);
410c58f1c22SToby Isaac 
4112cbb9b06SVaclav Hapla /*E
412af27ebaaSBarry Smith    DMCopyLabelsMode - Determines how `DMCopyLabels()` behaves when there is a `DMLabel` in the source and destination `DM`s with the same name
4132cbb9b06SVaclav Hapla 
41416a05f60SBarry Smith    Values:
41516a05f60SBarry Smith +  `DM_COPY_LABELS_REPLACE` - replace label in destination by label from source
41616a05f60SBarry Smith .  `DM_COPY_LABELS_KEEP`    - keep destination label
417af27ebaaSBarry Smith -  `DM_COPY_LABELS_FAIL`    - generate an error
41816a05f60SBarry Smith 
4192cbb9b06SVaclav Hapla    Level: advanced
4202cbb9b06SVaclav Hapla 
421af27ebaaSBarry Smith .seealso: [](ch_dmbase), `DMLabel`, `DM`, `DMCompareLabels()`, `DMRemoveLabel()`
4222cbb9b06SVaclav Hapla E*/
4239371c9d4SSatish Balay typedef enum {
4249371c9d4SSatish Balay   DM_COPY_LABELS_REPLACE,
4259371c9d4SSatish Balay   DM_COPY_LABELS_KEEP,
4269371c9d4SSatish Balay   DM_COPY_LABELS_FAIL
4279371c9d4SSatish Balay } DMCopyLabelsMode;
4282cbb9b06SVaclav Hapla PETSC_EXTERN const char *const DMCopyLabelsModes[];
4292cbb9b06SVaclav Hapla 
430c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetNumLabels(DM, PetscInt *);
431ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode DMGetLabelName(DM, PetscInt, const char *[]);
432c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMHasLabel(DM, const char[], PetscBool *);
433c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabel(DM, const char *, DMLabel *);
4344a7ee7d0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetLabel(DM, DMLabel);
435c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelByNum(DM, PetscInt, DMLabel *);
436c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMAddLabel(DM, DMLabel);
437c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMRemoveLabel(DM, const char[], DMLabel *);
438306894acSVaclav Hapla PETSC_EXTERN PetscErrorCode DMRemoveLabelBySelf(DM, DMLabel *, PetscBool);
439ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode DMCopyLabels(DM, DM, PetscCopyMode, PetscBool, DMCopyLabelsMode);
440ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode DMCompareLabels(DM, DM, PetscBool *, char *[]);
441c58f1c22SToby Isaac 
44257d50842SBarry Smith PETSC_EXTERN PetscErrorCode DMAddBoundary(DM, DMBoundaryConditionType, const char[], DMLabel, PetscInt, const PetscInt[], PetscInt, PetscInt, const PetscInt[], PetscVoidFn *, PetscVoidFn *, void *, PetscInt *);
443a6ba4734SToby Isaac PETSC_EXTERN PetscErrorCode DMIsBoundaryPoint(DM, PetscInt, PetscBool *);
44401468941SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasBound(DM, PetscBool *);
4454d6f44ffSToby Isaac 
4460709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMProjectFunction(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
4470709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMProjectFunctionLocal(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
4482c53366bSMatthew 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);
4491c531cf8SMatthew 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);
450191494d9SMatthew 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);
451ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode DMProjectFieldLabel(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);
4521c531cf8SMatthew 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);
453ece3a9fcSMatthew 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);
4540709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2Diff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
455b698f381SToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2GradientDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal[], PetscReal *);
4561189c1efSToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2FieldDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *);
4572e4af2aeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMComputeError(DM, Vec, PetscReal[], Vec *);
458ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasBasisTransform(DM, PetscBool *);
459ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyTransform(DM, DM);
4608320bc6fSPatrick Sanan 
4618320bc6fSPatrick Sanan PETSC_EXTERN PetscErrorCode DMGetCompatibility(DM, DM, PetscBool *, PetscBool *);
462c0f0dcc3SMatthew G. Knepley 
463*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMMonitorSet(DM, PetscErrorCode (*)(DM, PetscCtx), PetscCtx, PetscCtxDestroyFn *);
464c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitorCancel(DM);
465*2a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode DMMonitorSetFromOptions(DM, const char[], const char[], const char[], PetscErrorCode (*)(DM, PetscCtx), PetscErrorCode (*)(DM, PetscViewerAndFormat *), PetscBool *);
466c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitor(DM);
467c0f0dcc3SMatthew G. Knepley 
DMPolytopeTypeIsHybrid(DMPolytopeType ct)468c306944fSJed Brown static inline PetscBool DMPolytopeTypeIsHybrid(DMPolytopeType ct)
469c306944fSJed Brown {
470c306944fSJed Brown   switch (ct) {
471c306944fSJed Brown   case DM_POLYTOPE_POINT_PRISM_TENSOR:
472c306944fSJed Brown   case DM_POLYTOPE_SEG_PRISM_TENSOR:
473c306944fSJed Brown   case DM_POLYTOPE_TRI_PRISM_TENSOR:
474c306944fSJed Brown   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
475c306944fSJed Brown     return PETSC_TRUE;
476c306944fSJed Brown   default:
477c306944fSJed Brown     return PETSC_FALSE;
478c306944fSJed Brown   }
479c306944fSJed Brown }
480c306944fSJed Brown 
DMPolytopeTypeGetDim(DMPolytopeType ct)481d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetDim(DMPolytopeType ct)
482d71ae5a4SJacob Faibussowitsch {
483412e9a14SMatthew G. Knepley   switch (ct) {
484d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
485d71ae5a4SJacob Faibussowitsch     return 0;
486412e9a14SMatthew G. Knepley   case DM_POLYTOPE_SEGMENT:
487d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
488d71ae5a4SJacob Faibussowitsch     return 1;
489412e9a14SMatthew G. Knepley   case DM_POLYTOPE_TRIANGLE:
490412e9a14SMatthew G. Knepley   case DM_POLYTOPE_QUADRILATERAL:
491d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
492476787b7SMatthew G. Knepley   case DM_POLYTOPE_UNKNOWN_FACE:
493d71ae5a4SJacob Faibussowitsch     return 2;
494412e9a14SMatthew G. Knepley   case DM_POLYTOPE_TETRAHEDRON:
495412e9a14SMatthew G. Knepley   case DM_POLYTOPE_HEXAHEDRON:
496412e9a14SMatthew G. Knepley   case DM_POLYTOPE_TRI_PRISM:
497412e9a14SMatthew G. Knepley   case DM_POLYTOPE_TRI_PRISM_TENSOR:
498412e9a14SMatthew G. Knepley   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
499d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
500476787b7SMatthew G. Knepley   case DM_POLYTOPE_UNKNOWN_CELL:
501d71ae5a4SJacob Faibussowitsch     return 3;
502d71ae5a4SJacob Faibussowitsch   default:
503d71ae5a4SJacob Faibussowitsch     return -1;
504412e9a14SMatthew G. Knepley   }
505412e9a14SMatthew G. Knepley }
506412e9a14SMatthew G. Knepley 
DMPolytopeTypeGetConeSize(DMPolytopeType ct)507d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetConeSize(DMPolytopeType ct)
508d71ae5a4SJacob Faibussowitsch {
509412e9a14SMatthew G. Knepley   switch (ct) {
510d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
511d71ae5a4SJacob Faibussowitsch     return 0;
512d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
513d71ae5a4SJacob Faibussowitsch     return 2;
514d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
515d71ae5a4SJacob Faibussowitsch     return 2;
516d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
517d71ae5a4SJacob Faibussowitsch     return 3;
518d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
519d71ae5a4SJacob Faibussowitsch     return 4;
520d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
521d71ae5a4SJacob Faibussowitsch     return 4;
522d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
523d71ae5a4SJacob Faibussowitsch     return 4;
524d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
525d71ae5a4SJacob Faibussowitsch     return 6;
526d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
527d71ae5a4SJacob Faibussowitsch     return 5;
528d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
529d71ae5a4SJacob Faibussowitsch     return 5;
530d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
531d71ae5a4SJacob Faibussowitsch     return 6;
532d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
533d71ae5a4SJacob Faibussowitsch     return 5;
534d71ae5a4SJacob Faibussowitsch   default:
535d71ae5a4SJacob Faibussowitsch     return -1;
536412e9a14SMatthew G. Knepley   }
537412e9a14SMatthew G. Knepley }
538412e9a14SMatthew G. Knepley 
DMPolytopeTypeGetNumVertices(DMPolytopeType ct)539d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetNumVertices(DMPolytopeType ct)
540d71ae5a4SJacob Faibussowitsch {
541412e9a14SMatthew G. Knepley   switch (ct) {
542d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
543d71ae5a4SJacob Faibussowitsch     return 1;
544d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
545d71ae5a4SJacob Faibussowitsch     return 2;
546d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
547d71ae5a4SJacob Faibussowitsch     return 2;
548d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
549d71ae5a4SJacob Faibussowitsch     return 3;
550d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
551d71ae5a4SJacob Faibussowitsch     return 4;
552d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
553d71ae5a4SJacob Faibussowitsch     return 4;
554d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
555d71ae5a4SJacob Faibussowitsch     return 4;
556d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
557d71ae5a4SJacob Faibussowitsch     return 8;
558d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
559d71ae5a4SJacob Faibussowitsch     return 6;
560d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
561d71ae5a4SJacob Faibussowitsch     return 6;
562d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
563d71ae5a4SJacob Faibussowitsch     return 8;
564d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
565d71ae5a4SJacob Faibussowitsch     return 5;
566d71ae5a4SJacob Faibussowitsch   default:
567d71ae5a4SJacob Faibussowitsch     return -1;
568412e9a14SMatthew G. Knepley   }
569412e9a14SMatthew G. Knepley }
570412e9a14SMatthew G. Knepley 
DMPolytopeTypeSimpleShape(PetscInt dim,PetscBool simplex)571d71ae5a4SJacob Faibussowitsch static inline DMPolytopeType DMPolytopeTypeSimpleShape(PetscInt dim, PetscBool simplex)
572d71ae5a4SJacob Faibussowitsch {
5739371c9d4SSatish Balay   return dim == 0 ? DM_POLYTOPE_POINT : (dim == 1 ? DM_POLYTOPE_SEGMENT : (dim == 2 ? (simplex ? DM_POLYTOPE_TRIANGLE : DM_POLYTOPE_QUADRILATERAL) : (dim == 3 ? (simplex ? DM_POLYTOPE_TETRAHEDRON : DM_POLYTOPE_HEXAHEDRON) : DM_POLYTOPE_UNKNOWN)));
5749318fe57SMatthew G. Knepley }
5759318fe57SMatthew G. Knepley 
DMPolytopeTypeGetNumArrangements(DMPolytopeType ct)57685036b15SMatthew G. Knepley static inline PetscInt DMPolytopeTypeGetNumArrangements(DMPolytopeType ct)
577d71ae5a4SJacob Faibussowitsch {
578b5a892a1SMatthew G. Knepley   switch (ct) {
579d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
580d71ae5a4SJacob Faibussowitsch     return 1;
581d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
582d71ae5a4SJacob Faibussowitsch     return 2;
583d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
584d71ae5a4SJacob Faibussowitsch     return 2;
585d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
586d71ae5a4SJacob Faibussowitsch     return 6;
587d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
588d71ae5a4SJacob Faibussowitsch     return 8;
589d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
590d71ae5a4SJacob Faibussowitsch     return 4;
591d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
592d71ae5a4SJacob Faibussowitsch     return 24;
593d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
594d71ae5a4SJacob Faibussowitsch     return 48;
595d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
596d71ae5a4SJacob Faibussowitsch     return 12;
597d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
598d71ae5a4SJacob Faibussowitsch     return 12;
599d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
600d71ae5a4SJacob Faibussowitsch     return 16;
601d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
602d71ae5a4SJacob Faibussowitsch     return 8;
603d71ae5a4SJacob Faibussowitsch   default:
604d71ae5a4SJacob Faibussowitsch     return -1;
605b5a892a1SMatthew G. Knepley   }
606b5a892a1SMatthew G. Knepley }
607b5a892a1SMatthew G. Knepley 
608b5a892a1SMatthew G. Knepley /* An arrangement is a face order combined with an orientation for each face */
DMPolytopeTypeGetArrangement(DMPolytopeType ct,PetscInt o)60985036b15SMatthew G. Knepley static inline const PetscInt *DMPolytopeTypeGetArrangement(DMPolytopeType ct, PetscInt o)
610d71ae5a4SJacob Faibussowitsch {
611ef8b56bfSJed Brown   static const PetscInt pntArr[1 * 2] = {0, 0};
612b5a892a1SMatthew G. Knepley   /* a: swap */
6139371c9d4SSatish Balay   static const PetscInt segArr[2 * 2 * 2] = {1, 0, 0, 0, /* -1: a */
6149371c9d4SSatish Balay                                              0, 0, 1, 0,
6159371c9d4SSatish Balay                                              /*  0: e */};
616b5a892a1SMatthew G. Knepley   /* a: swap first two
617b5a892a1SMatthew G. Knepley      b: swap last two */
6189371c9d4SSatish Balay   static const PetscInt triArr[6 * 3 * 2] = {0, -1, 2, -1, 1, -1, /* -3: b */
619b5a892a1SMatthew G. Knepley                                              2, -1, 1, -1, 0, -1, /* -2: aba */
620b5a892a1SMatthew G. Knepley                                              1, -1, 0, -1, 2, -1, /* -1: a */
621b5a892a1SMatthew G. Knepley                                              0, 0,  1, 0,  2, 0,  /*  0: identity */
622b5a892a1SMatthew G. Knepley                                              1, 0,  2, 0,  0, 0,  /*  1: ba */
6239371c9d4SSatish Balay                                              2, 0,  0, 0,  1, 0,
6249371c9d4SSatish Balay                                              /*  2: ab */};
625b5a892a1SMatthew G. Knepley   /* a: forward cyclic permutation
626b5a892a1SMatthew G. Knepley      b: swap first and last pairs */
6279371c9d4SSatish Balay   static const PetscInt quadArr[8 * 4 * 2] = {1, -1, 0, -1, 3, -1, 2, -1, /* -4: b */
628b5a892a1SMatthew G. Knepley                                               0, -1, 3, -1, 2, -1, 1, -1, /* -3: b a^3 = a b */
629b5a892a1SMatthew G. Knepley                                               3, -1, 2, -1, 1, -1, 0, -1, /* -2: b a^2 = a^2 b */
630b5a892a1SMatthew G. Knepley                                               2, -1, 1, -1, 0, -1, 3, -1, /* -1: b a   = a^3 b */
631b5a892a1SMatthew G. Knepley                                               0, 0,  1, 0,  2, 0,  3, 0,  /*  0: identity */
632b5a892a1SMatthew G. Knepley                                               1, 0,  2, 0,  3, 0,  0, 0,  /*  1: a */
633b5a892a1SMatthew G. Knepley                                               2, 0,  3, 0,  0, 0,  1, 0,  /*  2: a^2 */
6349371c9d4SSatish Balay                                               3, 0,  0, 0,  1, 0,  2, 0,
6359371c9d4SSatish Balay                                               /*  3: a^3 */};
636b5a892a1SMatthew G. Knepley   /* r: rotate 180
637b5a892a1SMatthew G. Knepley      b: swap top and bottom segments */
6389371c9d4SSatish Balay   static const PetscInt tsegArr[4 * 4 * 2] = {1, -1, 0, -1, 3, -1, 2, -1, /* -2: r b */
639b5a892a1SMatthew G. Knepley                                               0, -1, 1, -1, 3, 0,  2, 0,  /* -1: r */
640b5a892a1SMatthew G. Knepley                                               0, 0,  1, 0,  2, 0,  3, 0,  /*  0: identity */
6419371c9d4SSatish Balay                                               1, 0,  0, 0,  2, -1, 3, -1,
6429371c9d4SSatish Balay                                               /*  1: b */};
643b5a892a1SMatthew G. Knepley   /* https://en.wikiversity.org/wiki/Symmetric_group_S4 */
6449371c9d4SSatish Balay   static const PetscInt tetArr[24 * 4 * 2] = {3, -2, 2, -3, 0, -1, 1, -1, /* -12: (1324)   p22 */
645b5a892a1SMatthew G. Knepley                                               3, -1, 1, -3, 2, -1, 0, -1, /* -11: (14)     p21 */
646b5a892a1SMatthew G. Knepley                                               3, -3, 0, -3, 1, -1, 2, -1, /* -10: (1234)   p18 */
647b5a892a1SMatthew G. Knepley                                               2, -1, 3, -1, 1, -3, 0, -2, /*  -9: (1423)   p17 */
648b5a892a1SMatthew G. Knepley                                               2, -3, 0, -1, 3, -2, 1, -3, /*  -8: (1342)   p13 */
649b5a892a1SMatthew G. Knepley                                               2, -2, 1, -2, 0, -2, 3, -2, /*  -7: (24)     p14 */
650b5a892a1SMatthew G. Knepley                                               1, -2, 0, -2, 2, -2, 3, -1, /*  -6: (34)     p6  */
651b5a892a1SMatthew G. Knepley                                               1, -1, 3, -3, 0, -3, 2, -2, /*  -5: (1243)   p10 */
652b5a892a1SMatthew G. Knepley                                               1, -3, 2, -1, 3, -1, 0, -3, /*  -4: (1432)   p9  */
653b5a892a1SMatthew G. Knepley                                               0, -3, 1, -1, 3, -3, 2, -3, /*  -3: (12)     p1  */
654b5a892a1SMatthew G. Knepley                                               0, -2, 2, -2, 1, -2, 3, -3, /*  -2: (23)     p2  */
655b5a892a1SMatthew G. Knepley                                               0, -1, 3, -2, 2, -3, 1, -2, /*  -1: (13)     p5  */
656b5a892a1SMatthew G. Knepley                                               0, 0,  1, 0,  2, 0,  3, 0,  /*   0: ()       p0  */
657b5a892a1SMatthew G. Knepley                                               0, 1,  3, 1,  1, 2,  2, 0,  /*   1: (123)    p4  */
658b5a892a1SMatthew G. Knepley                                               0, 2,  2, 1,  3, 0,  1, 2,  /*   2: (132)    p3  */
659b5a892a1SMatthew G. Knepley                                               1, 2,  0, 1,  3, 1,  2, 2,  /*   3: (12)(34) p7  */
660b5a892a1SMatthew G. Knepley                                               1, 0,  2, 0,  0, 0,  3, 1,  /*   4: (243)    p8  */
661b5a892a1SMatthew G. Knepley                                               1, 1,  3, 2,  2, 2,  0, 0,  /*   5: (143)    p11 */
662b5a892a1SMatthew G. Knepley                                               2, 1,  3, 0,  0, 2,  1, 0,  /*   6: (13)(24) p16 */
663b5a892a1SMatthew G. Knepley                                               2, 2,  1, 1,  3, 2,  0, 2,  /*   7: (142)    p15 */
664b5a892a1SMatthew G. Knepley                                               2, 0,  0, 0,  1, 0,  3, 2,  /*   8: (234)    p12 */
665b5a892a1SMatthew G. Knepley                                               3, 2,  2, 2,  1, 1,  0, 1,  /*   9: (14)(23) p23 */
666b5a892a1SMatthew G. Knepley                                               3, 0,  0, 2,  2, 1,  1, 1,  /*  10: (134)    p19 */
667b5a892a1SMatthew G. Knepley                                               3, 1,  1, 2,  0, 1,  2, 1 /*  11: (124)    p20 */};
668b5a892a1SMatthew G. Knepley   /* Each rotation determines a permutation of the four diagonals, and this defines the isomorphism with S_4 */
669ef8b56bfSJed Brown   static const PetscInt hexArr[48 * 6 * 2] = {
670b5a892a1SMatthew G. Knepley     2, -3, 3, -2, 4, -2, 5, -3, 1, -3, 0, -1, /* -24: reflect bottom and use -3 on top */
671b5a892a1SMatthew G. Knepley     4, -2, 5, -2, 0, -1, 1, -4, 3, -2, 2, -3, /* -23: reflect bottom and use -3 on top */
672b5a892a1SMatthew G. Knepley     5, -3, 4, -1, 1, -2, 0, -3, 3, -4, 2, -1, /* -22: reflect bottom and use -3 on top */
673b5a892a1SMatthew G. Knepley     3, -1, 2, -4, 4, -4, 5, -1, 0, -4, 1, -4, /* -21: reflect bottom and use -3 on top */
674b5a892a1SMatthew G. Knepley     3, -3, 2, -2, 5, -1, 4, -4, 1, -1, 0, -3, /* -20: reflect bottom and use -3 on top */
675b5a892a1SMatthew G. Knepley     4, -4, 5, -4, 1, -4, 0, -1, 2, -4, 3, -1, /* -19: reflect bottom and use -3 on top */
676b5a892a1SMatthew G. Knepley     2, -1, 3, -4, 5, -3, 4, -2, 0, -2, 1, -2, /* -18: reflect bottom and use -3 on top */
677b5a892a1SMatthew G. Knepley     5, -1, 4, -3, 0, -3, 1, -2, 2, -2, 3, -3, /* -17: reflect bottom and use -3 on top */
678b5a892a1SMatthew G. Knepley     4, -3, 5, -1, 3, -2, 2, -4, 1, -4, 0, -4, /* -16: reflect bottom and use -3 on top */
679b5a892a1SMatthew G. Knepley     5, -4, 4, -4, 3, -4, 2, -2, 0, -3, 1, -1, /* -15: reflect bottom and use -3 on top */
680b5a892a1SMatthew G. Knepley     3, -4, 2, -1, 1, -1, 0, -4, 4, -4, 5, -4, /* -14: reflect bottom and use -3 on top */
681b5a892a1SMatthew G. Knepley     2, -2, 3, -3, 0, -2, 1, -3, 4, -2, 5, -2, /* -13: reflect bottom and use -3 on top */
682b5a892a1SMatthew G. Knepley     1, -3, 0, -1, 4, -1, 5, -4, 3, -1, 2, -4, /* -12: reflect bottom and use -3 on top */
683b5a892a1SMatthew G. Knepley     1, -1, 0, -3, 5, -4, 4, -1, 2, -1, 3, -4, /* -11: reflect bottom and use -3 on top */
684b5a892a1SMatthew G. Knepley     5, -2, 4, -2, 2, -2, 3, -4, 1, -2, 0, -2, /* -10: reflect bottom and use -3 on top */
685b5a892a1SMatthew G. Knepley     1, -2, 0, -2, 2, -1, 3, -1, 4, -1, 5, -3, /*  -9: reflect bottom and use -3 on top */
686b5a892a1SMatthew G. Knepley     4, -1, 5, -3, 2, -4, 3, -2, 0, -1, 1, -3, /*  -8: reflect bottom and use -3 on top */
687b5a892a1SMatthew G. Knepley     3, -2, 2, -3, 0, -4, 1, -1, 5, -1, 4, -3, /*  -7: reflect bottom and use -3 on top */
688b5a892a1SMatthew G. Knepley     1, -4, 0, -4, 3, -1, 2, -1, 5, -4, 4, -4, /*  -6: reflect bottom and use -3 on top */
689b5a892a1SMatthew G. Knepley     2, -4, 3, -1, 1, -3, 0, -2, 5, -3, 4, -1, /*  -5: reflect bottom and use -3 on top */
690b5a892a1SMatthew G. Knepley     0, -4, 1, -4, 4, -3, 5, -2, 2, -3, 3, -2, /*  -4: reflect bottom and use -3 on top */
691b5a892a1SMatthew G. Knepley     0, -3, 1, -1, 3, -3, 2, -3, 4, -3, 5, -1, /*  -3: reflect bottom and use -3 on top */
692b5a892a1SMatthew G. Knepley     0, -2, 1, -2, 5, -2, 4, -3, 3, -3, 2, -2, /*  -2: reflect bottom and use -3 on top */
693b5a892a1SMatthew G. Knepley     0, -1, 1, -3, 2, -3, 3, -3, 5, -2, 4, -2, /*  -1: reflect bottom and use -3 on top */
694b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  5, 0,  /*   0: identity */
695b5a892a1SMatthew G. Knepley     0, 1,  1, 3,  5, 3,  4, 0,  2, 0,  3, 1,  /*   1: 90  rotation about z */
696b5a892a1SMatthew G. Knepley     0, 2,  1, 2,  3, 0,  2, 0,  5, 3,  4, 1,  /*   2: 180 rotation about z */
697b5a892a1SMatthew G. Knepley     0, 3,  1, 1,  4, 0,  5, 3,  3, 0,  2, 1,  /*   3: 270 rotation about z */
698b5a892a1SMatthew G. Knepley     2, 3,  3, 2,  1, 0,  0, 3,  4, 3,  5, 1,  /*   4: 90  rotation about x */
699b5a892a1SMatthew G. Knepley     1, 3,  0, 1,  3, 2,  2, 2,  4, 2,  5, 2,  /*   5: 180 rotation about x */
700b5a892a1SMatthew G. Knepley     3, 1,  2, 0,  0, 1,  1, 2,  4, 1,  5, 3,  /*   6: 270 rotation about x */
701b5a892a1SMatthew G. Knepley     4, 0,  5, 0,  2, 1,  3, 3,  1, 1,  0, 3,  /*   7: 90  rotation about y */
702b5a892a1SMatthew G. Knepley     1, 1,  0, 3,  2, 2,  3, 2,  5, 1,  4, 3,  /*   8: 180 rotation about y */
703b5a892a1SMatthew G. Knepley     5, 1,  4, 3,  2, 3,  3, 1,  0, 0,  1, 0,  /*   9: 270 rotation about y */
704b5a892a1SMatthew G. Knepley     1, 0,  0, 0,  5, 1,  4, 2,  3, 2,  2, 3,  /*  10: 180 rotation about x+y */
705b5a892a1SMatthew G. Knepley     1, 2,  0, 2,  4, 2,  5, 1,  2, 2,  3, 3,  /*  11: 180 rotation about x-y */
706b5a892a1SMatthew G. Knepley     2, 1,  3, 0,  0, 3,  1, 0,  5, 0,  4, 0,  /*  12: 180 rotation about y+z */
707b5a892a1SMatthew G. Knepley     3, 3,  2, 2,  1, 2,  0, 1,  5, 2,  4, 2,  /*  13: 180 rotation about y-z */
708b5a892a1SMatthew G. Knepley     5, 3,  4, 1,  3, 1,  2, 3,  1, 3,  0, 1,  /*  14: 180 rotation about z+x */
709b5a892a1SMatthew G. Knepley     4, 2,  5, 2,  3, 3,  2, 1,  0, 2,  1, 2,  /*  15: 180 rotation about z-x */
710b5a892a1SMatthew G. Knepley     5, 0,  4, 0,  0, 0,  1, 3,  3, 1,  2, 0,  /*  16: 120 rotation about x+y+z (v0v6) */
711b5a892a1SMatthew G. Knepley     2, 0,  3, 1,  5, 0,  4, 3,  1, 0,  0, 0,  /*  17: 240 rotation about x+y+z (v0v6) */
712b5a892a1SMatthew G. Knepley     4, 3,  5, 1,  1, 1,  0, 2,  3, 3,  2, 2,  /*  18: 120 rotation about x+y-z (v4v2) */
713b5a892a1SMatthew G. Knepley     3, 2,  2, 3,  5, 2,  4, 1,  0, 1,  1, 3,  /*  19: 240 rotation about x+y-z (v4v2) */
714b5a892a1SMatthew G. Knepley     3, 0,  2, 1,  4, 1,  5, 2,  1, 2,  0, 2,  /*  20: 120 rotation about x-y+z (v1v5) */
715b5a892a1SMatthew G. Knepley     5, 2,  4, 2,  1, 3,  0, 0,  2, 3,  3, 2,  /*  21: 240 rotation about x-y+z (v1v5) */
716b5a892a1SMatthew G. Knepley     4, 1,  5, 3,  0, 2,  1, 1,  2, 1,  3, 0,  /*  22: 120 rotation about x-y-z (v7v3) */
717b5a892a1SMatthew G. Knepley     2, 2,  3, 3,  4, 3,  5, 0,  0, 3,  1, 1,  /*  23: 240 rotation about x-y-z (v7v3) */
718b5a892a1SMatthew G. Knepley   };
719ef8b56bfSJed Brown   static const PetscInt tripArr[12 * 5 * 2] = {
720b5a892a1SMatthew G. Knepley     1, -3, 0, -1, 3, -1, 4, -1, 2, -1, /* -6: reflect bottom and top */
721b5a892a1SMatthew G. Knepley     1, -1, 0, -3, 4, -1, 2, -1, 3, -1, /* -5: reflect bottom and top */
722b5a892a1SMatthew G. Knepley     1, -2, 0, -2, 2, -1, 3, -1, 4, -1, /* -4: reflect bottom and top */
723b5a892a1SMatthew G. Knepley     0, -3, 1, -1, 3, -3, 2, -3, 4, -3, /* -3: reflect bottom and top */
724b5a892a1SMatthew G. Knepley     0, -2, 1, -2, 4, -3, 3, -3, 2, -3, /* -2: reflect bottom and top */
725b5a892a1SMatthew G. Knepley     0, -1, 1, -3, 2, -3, 4, -3, 3, -3, /* -1: reflect bottom and top */
726b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  /*  0: identity */
727b5a892a1SMatthew G. Knepley     0, 1,  1, 2,  4, 0,  2, 0,  3, 0,  /*  1: 120 rotation about z */
728b5a892a1SMatthew G. Knepley     0, 2,  1, 1,  3, 0,  4, 0,  2, 0,  /*  2: 240 rotation about z */
729b5a892a1SMatthew G. Knepley     1, 1,  0, 2,  2, 2,  4, 2,  3, 2,  /*  3: 180 rotation about y of 0 */
730b5a892a1SMatthew G. Knepley     1, 0,  0, 0,  4, 2,  3, 2,  2, 2,  /*  4: 180 rotation about y of 1 */
731b5a892a1SMatthew G. Knepley     1, 2,  0, 1,  3, 2,  2, 2,  4, 2,  /*  5: 180 rotation about y of 2 */
732b5a892a1SMatthew G. Knepley   };
733b5a892a1SMatthew G. Knepley   /* a: rotate 120 about z
734b5a892a1SMatthew G. Knepley      b: swap top and bottom segments
735b5a892a1SMatthew G. Knepley      r: reflect */
736ef8b56bfSJed Brown   static const PetscInt ttriArr[12 * 5 * 2] = {
737b5a892a1SMatthew G. Knepley     1, -3, 0, -3, 2, -2, 4, -2, 3, -2, /* -6: r b a^2 */
738b5a892a1SMatthew G. Knepley     1, -2, 0, -2, 4, -2, 3, -2, 2, -2, /* -5: r b a */
739b5a892a1SMatthew G. Knepley     1, -1, 0, -1, 3, -2, 2, -2, 4, -2, /* -4: r b */
740b5a892a1SMatthew G. Knepley     0, -3, 1, -3, 2, -1, 4, -1, 3, -1, /* -3: r a^2 */
741b5a892a1SMatthew G. Knepley     0, -2, 1, -2, 4, -1, 3, -1, 2, -1, /* -2: r a */
742b5a892a1SMatthew G. Knepley     0, -1, 1, -1, 3, -1, 2, -1, 4, -1, /* -1: r */
743b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  /*  0: identity */
744b5a892a1SMatthew G. Knepley     0, 1,  1, 1,  3, 0,  4, 0,  2, 0,  /*  1: a */
745b5a892a1SMatthew G. Knepley     0, 2,  1, 2,  4, 0,  2, 0,  3, 0,  /*  2: a^2 */
746b5a892a1SMatthew G. Knepley     1, 0,  0, 0,  2, 1,  3, 1,  4, 1,  /*  3: b */
747b5a892a1SMatthew G. Knepley     1, 1,  0, 1,  3, 1,  4, 1,  2, 1,  /*  4: b a */
748b5a892a1SMatthew G. Knepley     1, 2,  0, 2,  4, 1,  2, 1,  3, 1,  /*  5: b a^2 */
749b5a892a1SMatthew G. Knepley   };
750b5a892a1SMatthew G. Knepley   /* a: rotate 90 about z
751b5a892a1SMatthew G. Knepley      b: swap top and bottom segments
752b5a892a1SMatthew G. Knepley      r: reflect */
753ef8b56bfSJed Brown   static const PetscInt tquadArr[16 * 6 * 2] = {
754b5a892a1SMatthew G. Knepley     1, -4, 0, -4, 3, -2, 2, -2, 5, -2, 4, -2, /* -8: r b a^3 */
755b5a892a1SMatthew G. Knepley     1, -3, 0, -3, 2, -2, 5, -2, 4, -2, 3, -2, /* -7: r b a^2 */
756b5a892a1SMatthew G. Knepley     1, -2, 0, -2, 5, -2, 4, -2, 3, -2, 2, -2, /* -6: r b a */
757b5a892a1SMatthew G. Knepley     1, -1, 0, -1, 4, -2, 3, -2, 2, -2, 5, -2, /* -5: r b */
758b5a892a1SMatthew G. Knepley     0, -4, 1, -4, 3, -1, 2, -1, 5, -1, 4, -1, /* -4: r a^3 */
759b5a892a1SMatthew G. Knepley     0, -3, 1, -3, 2, -1, 5, -1, 4, -1, 3, -1, /* -3: r a^2 */
760b5a892a1SMatthew G. Knepley     0, -2, 1, -2, 5, -1, 4, -1, 3, -1, 2, -1, /* -2: r a */
761b5a892a1SMatthew G. Knepley     0, -1, 1, -1, 4, -1, 3, -1, 2, -1, 5, -1, /* -1: r */
762b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  5, 0,  /*  0: identity */
763b5a892a1SMatthew G. Knepley     0, 1,  1, 1,  3, 0,  4, 0,  5, 0,  2, 0,  /*  1: a */
764b5a892a1SMatthew G. Knepley     0, 2,  1, 2,  4, 0,  5, 0,  2, 0,  3, 0,  /*  2: a^2 */
765b5a892a1SMatthew G. Knepley     0, 3,  1, 3,  5, 0,  2, 0,  3, 0,  4, 0,  /*  3: a^3 */
766b5a892a1SMatthew G. Knepley     1, 0,  0, 0,  2, 1,  3, 1,  4, 1,  5, 1,  /*  4: b */
767b5a892a1SMatthew G. Knepley     1, 1,  0, 1,  3, 1,  4, 1,  5, 1,  2, 1,  /*  5: b a */
768b5a892a1SMatthew G. Knepley     1, 2,  0, 2,  4, 1,  5, 1,  2, 1,  3, 1,  /*  6: b a^2 */
769b5a892a1SMatthew G. Knepley     1, 3,  0, 3,  5, 1,  2, 1,  3, 1,  4, 1,  /*  7: b a^3 */
770b5a892a1SMatthew G. Knepley   };
771ef8b56bfSJed Brown   static const PetscInt pyrArr[8 * 5 * 2] = {
772b5a892a1SMatthew G. Knepley     0, -4, 2, -3, 1, -3, 4, -3, 3, -3, /* -4: Reflect bottom face */
773b5a892a1SMatthew G. Knepley     0, -3, 3, -3, 2, -3, 1, -3, 4, -3, /* -3: Reflect bottom face */
774b5a892a1SMatthew G. Knepley     0, -2, 4, -3, 3, -3, 2, -3, 1, -3, /* -2: Reflect bottom face */
775b5a892a1SMatthew G. Knepley     0, -1, 1, -3, 4, -3, 3, -3, 2, -3, /* -1: Reflect bottom face */
776b5a892a1SMatthew G. Knepley     0, 0,  1, 0,  2, 0,  3, 0,  4, 0,  /*  0: identity */
777b5a892a1SMatthew G. Knepley     0, 1,  4, 0,  1, 0,  2, 0,  3, 0,  /*  1:  90 rotation about z */
778b5a892a1SMatthew G. Knepley     0, 2,  3, 0,  4, 0,  1, 0,  2, 0,  /*  2: 180 rotation about z */
779b5a892a1SMatthew G. Knepley     0, 3,  2, 0,  3, 0,  4, 0,  1, 0,  /*  3: 270 rotation about z */
780b5a892a1SMatthew G. Knepley   };
781b5a892a1SMatthew G. Knepley   switch (ct) {
782d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
783d71ae5a4SJacob Faibussowitsch     return pntArr;
784d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
785d71ae5a4SJacob Faibussowitsch     return &segArr[(o + 1) * 2 * 2];
786d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
787d71ae5a4SJacob Faibussowitsch     return &segArr[(o + 1) * 2 * 2];
788d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
789d71ae5a4SJacob Faibussowitsch     return &triArr[(o + 3) * 3 * 2];
790d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
791d71ae5a4SJacob Faibussowitsch     return &quadArr[(o + 4) * 4 * 2];
792d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
793d71ae5a4SJacob Faibussowitsch     return &tsegArr[(o + 2) * 4 * 2];
794d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
795d71ae5a4SJacob Faibussowitsch     return &tetArr[(o + 12) * 4 * 2];
796d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
797d71ae5a4SJacob Faibussowitsch     return &hexArr[(o + 24) * 6 * 2];
798d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
799d71ae5a4SJacob Faibussowitsch     return &tripArr[(o + 6) * 5 * 2];
800d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
801d71ae5a4SJacob Faibussowitsch     return &ttriArr[(o + 6) * 5 * 2];
802d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
803d71ae5a4SJacob Faibussowitsch     return &tquadArr[(o + 8) * 6 * 2];
804d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
805d71ae5a4SJacob Faibussowitsch     return &pyrArr[(o + 4) * 5 * 2];
806d71ae5a4SJacob Faibussowitsch   default:
807f22e26b7SPierre Jolivet     return PETSC_NULLPTR;
808b5a892a1SMatthew G. Knepley   }
809b5a892a1SMatthew G. Knepley }
810b5a892a1SMatthew G. Knepley 
811d5b43468SJose E. Roman /* A vertex arrangement is a vertex order */
DMPolytopeTypeGetVertexArrangement(DMPolytopeType ct,PetscInt o)81285036b15SMatthew G. Knepley static inline const PetscInt *DMPolytopeTypeGetVertexArrangement(DMPolytopeType ct, PetscInt o)
813d71ae5a4SJacob Faibussowitsch {
814ef8b56bfSJed Brown   static const PetscInt pntVerts[1]      = {0};
8159371c9d4SSatish Balay   static const PetscInt segVerts[2 * 2]  = {1, 0, 0, 1};
8169371c9d4SSatish Balay   static const PetscInt triVerts[6 * 3]  = {1, 0, 2, 0, 2, 1, 2, 1, 0, 0, 1, 2, 1, 2, 0, 2, 0, 1};
8179371c9d4SSatish Balay   static const PetscInt quadVerts[8 * 4] = {2, 1, 0, 3, 1, 0, 3, 2, 0, 3, 2, 1, 3, 2, 1, 0, 0, 1, 2, 3, 1, 2, 3, 0, 2, 3, 0, 1, 3, 0, 1, 2};
8189371c9d4SSatish Balay   static const PetscInt tsegVerts[4 * 4] = {3, 2, 1, 0, 1, 0, 3, 2, 0, 1, 2, 3, 2, 3, 0, 1};
8199371c9d4SSatish Balay   static const PetscInt tetVerts[24 * 4] = {2, 3, 1, 0, /* -12: (1324)   p22 */
820b5a892a1SMatthew G. Knepley                                             3, 1, 2, 0, /* -11: (14)     p21 */
821b5a892a1SMatthew G. Knepley                                             1, 2, 3, 0, /* -10: (1234)   p18 */
822b5a892a1SMatthew G. Knepley                                             3, 2, 0, 1, /*  -9: (1423)   p17 */
823b5a892a1SMatthew G. Knepley                                             2, 0, 3, 1, /*  -8: (1342)   p13 */
824b5a892a1SMatthew G. Knepley                                             0, 3, 2, 1, /*  -7: (24)     p14 */
825b5a892a1SMatthew G. Knepley                                             0, 1, 3, 2, /*  -6: (34)     p6  */
826b5a892a1SMatthew G. Knepley                                             1, 3, 0, 2, /*  -5: (1243)   p10 */
827b5a892a1SMatthew G. Knepley                                             3, 0, 1, 2, /*  -4: (1432    p9  */
828b5a892a1SMatthew G. Knepley                                             1, 0, 2, 3, /*  -3: (12)     p1  */
829b5a892a1SMatthew G. Knepley                                             0, 2, 1, 3, /*  -2: (23)     p2  */
830b5a892a1SMatthew G. Knepley                                             2, 1, 0, 3, /*  -1: (13)     p5  */
831b5a892a1SMatthew G. Knepley                                             0, 1, 2, 3, /*   0: ()       p0  */
832b5a892a1SMatthew G. Knepley                                             1, 2, 0, 3, /*   1: (123)    p4  */
833b5a892a1SMatthew G. Knepley                                             2, 0, 1, 3, /*   2: (132)    p3  */
834b5a892a1SMatthew G. Knepley                                             1, 0, 3, 2, /*   3: (12)(34) p7  */
835b5a892a1SMatthew G. Knepley                                             0, 3, 1, 2, /*   4: (243)    p8  */
836b5a892a1SMatthew G. Knepley                                             3, 1, 0, 2, /*   5: (143)    p11 */
837b5a892a1SMatthew G. Knepley                                             2, 3, 0, 1, /*   6: (13)(24) p16 */
838b5a892a1SMatthew G. Knepley                                             3, 0, 2, 1, /*   7: (142)    p15 */
839b5a892a1SMatthew G. Knepley                                             0, 2, 3, 1, /*   8: (234)    p12 */
840b5a892a1SMatthew G. Knepley                                             3, 2, 1, 0, /*   9: (14)(23) p23 */
841b5a892a1SMatthew G. Knepley                                             2, 1, 3, 0, /*  10: (134)    p19 */
842b5a892a1SMatthew G. Knepley                                             1, 3, 2, 0 /*  11: (124)    p20 */};
843ef8b56bfSJed Brown   static const PetscInt hexVerts[48 * 8] = {
844b5a892a1SMatthew G. Knepley     3, 0, 4, 5, 2, 6, 7, 1, /* -24: reflected 23 */
845b5a892a1SMatthew G. Knepley     3, 5, 6, 2, 0, 1, 7, 4, /* -23: reflected 22 */
846b5a892a1SMatthew G. Knepley     4, 0, 1, 7, 5, 6, 2, 3, /* -22: reflected 21 */
847b5a892a1SMatthew G. Knepley     6, 7, 1, 2, 5, 3, 0, 4, /* -21: reflected 20 */
848b5a892a1SMatthew G. Knepley     1, 2, 6, 7, 0, 4, 5, 3, /* -20: reflected 19 */
849b5a892a1SMatthew G. Knepley     6, 2, 3, 5, 7, 4, 0, 1, /* -19: reflected 18 */
850b5a892a1SMatthew G. Knepley     4, 5, 3, 0, 7, 1, 2, 6, /* -18: reflected 17 */
851b5a892a1SMatthew G. Knepley     1, 7, 4, 0, 2, 3, 5, 6, /* -17: reflected 16 */
852b5a892a1SMatthew G. Knepley     2, 3, 5, 6, 1, 7, 4, 0, /* -16: reflected 15 */
853b5a892a1SMatthew G. Knepley     7, 4, 0, 1, 6, 2, 3, 5, /* -15: reflected 14 */
854b5a892a1SMatthew G. Knepley     7, 1, 2, 6, 4, 5, 3, 0, /* -14: reflected 13 */
855b5a892a1SMatthew G. Knepley     0, 4, 5, 3, 1, 2, 6, 7, /* -13: reflected 12 */
856b5a892a1SMatthew G. Knepley     5, 4, 7, 6, 3, 2, 1, 0, /* -12: reflected 11 */
857b5a892a1SMatthew G. Knepley     7, 6, 5, 4, 1, 0, 3, 2, /* -11: reflected 10 */
858b5a892a1SMatthew G. Knepley     0, 1, 7, 4, 3, 5, 6, 2, /* -10: reflected  9 */
859b5a892a1SMatthew G. Knepley     4, 7, 6, 5, 0, 3, 2, 1, /*  -9: reflected  8 */
860b5a892a1SMatthew G. Knepley     5, 6, 2, 3, 4, 0, 1, 7, /*  -8: reflected  7 */
861b5a892a1SMatthew G. Knepley     2, 6, 7, 1, 3, 0, 4, 5, /*  -7: reflected  6 */
862b5a892a1SMatthew G. Knepley     6, 5, 4, 7, 2, 1, 0, 3, /*  -6: reflected  5 */
863b5a892a1SMatthew G. Knepley     5, 3, 0, 4, 6, 7, 1, 2, /*  -5: reflected  4 */
864b5a892a1SMatthew G. Knepley     2, 1, 0, 3, 6, 5, 4, 7, /*  -4: reflected  3 */
865b5a892a1SMatthew G. Knepley     1, 0, 3, 2, 7, 6, 5, 4, /*  -3: reflected  2 */
866b5a892a1SMatthew G. Knepley     0, 3, 2, 1, 4, 7, 6, 5, /*  -2: reflected  1 */
867b5a892a1SMatthew G. Knepley     3, 2, 1, 0, 5, 4, 7, 6, /*  -1: reflected  0 */
868b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, 5, 6, 7, /*   0: identity */
869b5a892a1SMatthew G. Knepley     1, 2, 3, 0, 7, 4, 5, 6, /*   1: 90  rotation about z */
870b5a892a1SMatthew G. Knepley     2, 3, 0, 1, 6, 7, 4, 5, /*   2: 180 rotation about z */
871b5a892a1SMatthew G. Knepley     3, 0, 1, 2, 5, 6, 7, 4, /*   3: 270 rotation about z */
872b5a892a1SMatthew G. Knepley     4, 0, 3, 5, 7, 6, 2, 1, /*   4: 90  rotation about x */
873b5a892a1SMatthew G. Knepley     7, 4, 5, 6, 1, 2, 3, 0, /*   5: 180 rotation about x */
874b5a892a1SMatthew G. Knepley     1, 7, 6, 2, 0, 3, 5, 4, /*   6: 270 rotation about x */
875b5a892a1SMatthew G. Knepley     3, 2, 6, 5, 0, 4, 7, 1, /*   7: 90  rotation about y */
876b5a892a1SMatthew G. Knepley     5, 6, 7, 4, 3, 0, 1, 2, /*   8: 180 rotation about y */
877b5a892a1SMatthew G. Knepley     4, 7, 1, 0, 5, 3, 2, 6, /*   9: 270 rotation about y */
878b5a892a1SMatthew G. Knepley     4, 5, 6, 7, 0, 1, 2, 3, /*  10: 180 rotation about x+y */
879b5a892a1SMatthew G. Knepley     6, 7, 4, 5, 2, 3, 0, 1, /*  11: 180 rotation about x-y */
880b5a892a1SMatthew G. Knepley     3, 5, 4, 0, 2, 1, 7, 6, /*  12: 180 rotation about y+z */
881b5a892a1SMatthew G. Knepley     6, 2, 1, 7, 5, 4, 0, 3, /*  13: 180 rotation about y-z */
882b5a892a1SMatthew G. Knepley     1, 0, 4, 7, 2, 6, 5, 3, /*  14: 180 rotation about z+x */
883b5a892a1SMatthew G. Knepley     6, 5, 3, 2, 7, 1, 0, 4, /*  15: 180 rotation about z-x */
884b5a892a1SMatthew G. Knepley     0, 4, 7, 1, 3, 2, 6, 5, /*  16: 120 rotation about x+y+z (v0v6) */
885b5a892a1SMatthew G. Knepley     0, 3, 5, 4, 1, 7, 6, 2, /*  17: 240 rotation about x+y+z (v0v6) */
886b5a892a1SMatthew G. Knepley     5, 3, 2, 6, 4, 7, 1, 0, /*  18: 120 rotation about x+y-z (v4v2) */
887b5a892a1SMatthew G. Knepley     7, 6, 2, 1, 4, 0, 3, 5, /*  19: 240 rotation about x+y-z (v4v2) */
888b5a892a1SMatthew G. Knepley     2, 1, 7, 6, 3, 5, 4, 0, /*  20: 120 rotation about x-y+z (v1v5) */
889b5a892a1SMatthew G. Knepley     7, 1, 0, 4, 6, 5, 3, 2, /*  21: 240 rotation about x-y+z (v1v5) */
890b5a892a1SMatthew G. Knepley     2, 6, 5, 3, 1, 0, 4, 7, /*  22: 120 rotation about x-y-z (v7v3) */
891b5a892a1SMatthew G. Knepley     5, 4, 0, 3, 6, 2, 1, 7, /*  23: 240 rotation about x-y-z (v7v3) */
892b5a892a1SMatthew G. Knepley   };
893ef8b56bfSJed Brown   static const PetscInt tripVerts[12 * 6] = {
894b5a892a1SMatthew G. Knepley     4, 3, 5, 2, 1, 0, /* -6: reflect bottom and top */
895b5a892a1SMatthew G. Knepley     5, 4, 3, 1, 0, 2, /* -5: reflect bottom and top */
896b5a892a1SMatthew G. Knepley     3, 5, 4, 0, 2, 1, /* -4: reflect bottom and top */
897b5a892a1SMatthew G. Knepley     1, 0, 2, 5, 4, 3, /* -3: reflect bottom and top */
898b5a892a1SMatthew G. Knepley     0, 2, 1, 3, 5, 4, /* -2: reflect bottom and top */
899b5a892a1SMatthew G. Knepley     2, 1, 0, 4, 3, 5, /* -1: reflect bottom and top */
900b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, 5, /*  0: identity */
901b5a892a1SMatthew G. Knepley     1, 2, 0, 5, 3, 4, /*  1: 120 rotation about z */
902b5a892a1SMatthew G. Knepley     2, 0, 1, 4, 5, 3, /*  2: 240 rotation about z */
903b5a892a1SMatthew G. Knepley     4, 5, 3, 2, 0, 1, /*  3: 180 rotation about y of 0 */
904b5a892a1SMatthew G. Knepley     3, 4, 5, 0, 1, 2, /*  4: 180 rotation about y of 1 */
905b5a892a1SMatthew G. Knepley     5, 3, 4, 1, 2, 0, /*  5: 180 rotation about y of 2 */
906b5a892a1SMatthew G. Knepley   };
907ef8b56bfSJed Brown   static const PetscInt ttriVerts[12 * 6] = {
908b5a892a1SMatthew G. Knepley     4, 3, 5, 1, 0, 2, /* -6: r b a^2 */
909b5a892a1SMatthew G. Knepley     3, 5, 4, 0, 2, 1, /* -5: r b a */
910b5a892a1SMatthew G. Knepley     5, 4, 3, 2, 1, 0, /* -4: r b */
911b5a892a1SMatthew G. Knepley     1, 0, 2, 4, 3, 5, /* -3: r a^2 */
912b5a892a1SMatthew G. Knepley     0, 2, 1, 3, 5, 4, /* -2: r a */
913b5a892a1SMatthew G. Knepley     2, 1, 0, 5, 4, 3, /* -1: r */
914b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, 5, /*  0: identity */
915b5a892a1SMatthew G. Knepley     1, 2, 0, 4, 5, 3, /*  1: a */
916b5a892a1SMatthew G. Knepley     2, 0, 1, 5, 3, 4, /*  2: a^2 */
917b5a892a1SMatthew G. Knepley     3, 4, 5, 0, 1, 2, /*  3: b */
918b5a892a1SMatthew G. Knepley     4, 5, 3, 1, 2, 0, /*  4: b a */
919b5a892a1SMatthew G. Knepley     5, 3, 4, 2, 0, 1, /*  5: b a^2 */
920b5a892a1SMatthew G. Knepley   };
921b5a892a1SMatthew G. Knepley   /* a: rotate 90 about z
922b5a892a1SMatthew G. Knepley      b: swap top and bottom segments
923b5a892a1SMatthew G. Knepley      r: reflect */
924ef8b56bfSJed Brown   static const PetscInt tquadVerts[16 * 8] = {
925b5a892a1SMatthew G. Knepley     6, 5, 4, 7, 2, 1, 0, 3, /* -8: r b a^3 */
926b5a892a1SMatthew G. Knepley     5, 4, 7, 6, 1, 0, 3, 2, /* -7: r b a^2 */
927b5a892a1SMatthew G. Knepley     4, 7, 6, 5, 0, 3, 2, 1, /* -6: r b a */
928b5a892a1SMatthew G. Knepley     7, 6, 5, 4, 3, 2, 1, 0, /* -5: r b */
929b5a892a1SMatthew G. Knepley     2, 1, 0, 3, 6, 5, 4, 7, /* -4: r a^3 */
930b5a892a1SMatthew G. Knepley     1, 0, 3, 2, 5, 4, 7, 6, /* -3: r a^2 */
931b5a892a1SMatthew G. Knepley     0, 3, 2, 1, 4, 7, 6, 5, /* -2: r a */
932b5a892a1SMatthew G. Knepley     3, 2, 1, 0, 7, 6, 5, 4, /* -1: r */
933b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, 5, 6, 7, /*  0: identity */
934b5a892a1SMatthew G. Knepley     1, 2, 3, 0, 5, 6, 7, 4, /*  1: a */
935b5a892a1SMatthew G. Knepley     2, 3, 0, 1, 6, 7, 4, 5, /*  2: a^2 */
936b5a892a1SMatthew G. Knepley     3, 0, 1, 2, 7, 4, 5, 6, /*  3: a^3 */
937b5a892a1SMatthew G. Knepley     4, 5, 6, 7, 0, 1, 2, 3, /*  4: b */
938b5a892a1SMatthew G. Knepley     5, 6, 7, 4, 1, 2, 3, 0, /*  5: b a */
939b5a892a1SMatthew G. Knepley     6, 7, 4, 5, 2, 3, 0, 1, /*  6: b a^2 */
940b5a892a1SMatthew G. Knepley     7, 4, 5, 6, 3, 0, 1, 2, /*  7: b a^3 */
941b5a892a1SMatthew G. Knepley   };
942ef8b56bfSJed Brown   static const PetscInt pyrVerts[8 * 5] = {
943b5a892a1SMatthew G. Knepley     2, 1, 0, 3, 4, /* -4: Reflect bottom face */
944b5a892a1SMatthew G. Knepley     1, 0, 3, 2, 4, /* -3: Reflect bottom face */
945b5a892a1SMatthew G. Knepley     0, 3, 2, 1, 4, /* -2: Reflect bottom face */
946b5a892a1SMatthew G. Knepley     3, 2, 1, 0, 4, /* -1: Reflect bottom face */
947b5a892a1SMatthew G. Knepley     0, 1, 2, 3, 4, /*  0: identity */
948b5a892a1SMatthew G. Knepley     1, 2, 3, 0, 4, /*  1:  90 rotation about z */
949b5a892a1SMatthew G. Knepley     2, 3, 0, 1, 4, /*  2: 180 rotation about z */
950b5a892a1SMatthew G. Knepley     3, 0, 1, 2, 4, /*  3: 270 rotation about z */
951b5a892a1SMatthew G. Knepley   };
952b5a892a1SMatthew G. Knepley   switch (ct) {
953d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
954d71ae5a4SJacob Faibussowitsch     return pntVerts;
955d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEGMENT:
956d71ae5a4SJacob Faibussowitsch     return &segVerts[(o + 1) * 2];
957d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
958d71ae5a4SJacob Faibussowitsch     return &segVerts[(o + 1) * 2];
959d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
960d71ae5a4SJacob Faibussowitsch     return &triVerts[(o + 3) * 3];
961d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
962d71ae5a4SJacob Faibussowitsch     return &quadVerts[(o + 4) * 4];
963d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
964d71ae5a4SJacob Faibussowitsch     return &tsegVerts[(o + 2) * 4];
965d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
966d71ae5a4SJacob Faibussowitsch     return &tetVerts[(o + 12) * 4];
967d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
968d71ae5a4SJacob Faibussowitsch     return &hexVerts[(o + 24) * 8];
969d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
970d71ae5a4SJacob Faibussowitsch     return &tripVerts[(o + 6) * 6];
971d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
972d71ae5a4SJacob Faibussowitsch     return &ttriVerts[(o + 6) * 6];
973d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
974d71ae5a4SJacob Faibussowitsch     return &tquadVerts[(o + 8) * 8];
975d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
976d71ae5a4SJacob Faibussowitsch     return &pyrVerts[(o + 4) * 5];
977d71ae5a4SJacob Faibussowitsch   default:
978f22e26b7SPierre Jolivet     return PETSC_NULLPTR;
979b5a892a1SMatthew G. Knepley   }
980b5a892a1SMatthew G. Knepley }
981b5a892a1SMatthew G. Knepley 
982b5a892a1SMatthew G. Knepley /* This is orientation o1 acting on orientation o2 */
DMPolytopeTypeComposeOrientation(DMPolytopeType ct,PetscInt o1,PetscInt o2)983d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeComposeOrientation(DMPolytopeType ct, PetscInt o1, PetscInt o2)
984d71ae5a4SJacob Faibussowitsch {
9859371c9d4SSatish Balay   static const PetscInt segMult[2 * 2]   = {0, -1, -1, 0};
9869371c9d4SSatish Balay   static const PetscInt triMult[6 * 6]   = {0, 2, 1, -3, -1, -2, 1, 0, 2, -2, -3, -1, 2, 1, 0, -1, -2, -3, -3, -2, -1, 0, 1, 2, -2, -1, -3, 1, 2, 0, -1, -3, -2, 2, 0, 1};
9879371c9d4SSatish Balay   static const PetscInt quadMult[8 * 8]  = {0,  3,  2,  1,  -4, -1, -2, -3, 1,  0,  3,  2,  -3, -4, -1, -2, 2,  1,  0,  3,  -2, -3, -4, -1, 3,  2,  1,  0,  -1, -2, -3, -4,
9889371c9d4SSatish Balay                                             -4, -3, -2, -1, 0,  1,  2,  3,  -3, -2, -1, -4, 1,  2,  3,  0,  -2, -1, -4, -3, 2,  3,  0,  1,  -1, -4, -3, -2, 3,  0,  1,  2};
9899371c9d4SSatish Balay   static const PetscInt tsegMult[4 * 4]  = {0, 1, -2, -1, 1, 0, -1, -2, -2, -1, 0, 1, -1, -2, 1, 0};
990ef8b56bfSJed Brown   static const PetscInt tetMult[24 * 24] = {
9919371c9d4SSatish Balay     3,   2,   7,   0,   5,   10,  9,   8,   1,   6,   11,  4,   -12, -7,  -5,  -9,  -10, -2,  -6,  -1,  -11, -3,  -4,  -8,  4,   0,   8,   1,   3,   11,  10,  6,   2,   7,   9,   5,   -11, -9,  -4,  -8,  -12, -1,  -5,  -3,  -10, -2,  -6,  -7,
9929371c9d4SSatish Balay     5,   1,   6,   2,   4,   9,   11,  7,   0,   8,   10,  3,   -10, -8,  -6,  -7,  -11, -3,  -4,  -2,  -12, -1,  -5,  -9,  0,   8,   4,   3,   11,  1,   6,   2,   10,  9,   5,   7,   -9,  -4,  -11, -12, -1,  -8,  -3,  -10, -5,  -6,  -7,  -2,
9939371c9d4SSatish Balay     1,   6,   5,   4,   9,   2,   7,   0,   11,  10,  3,   8,   -8,  -6,  -10, -11, -3,  -7,  -2,  -12, -4,  -5,  -9,  -1,  2,   7,   3,   5,   10,  0,   8,   1,   9,   11,  4,   6,   -7,  -5,  -12, -10, -2,  -9,  -1,  -11, -6,  -4,  -8,  -3,
9949371c9d4SSatish Balay     6,   5,   1,   9,   2,   4,   0,   11,  7,   3,   8,   10,  -6,  -10, -8,  -3,  -7,  -11, -12, -4,  -2,  -9,  -1,  -5,  7,   3,   2,   10,  0,   5,   1,   9,   8,   4,   6,   11,  -5,  -12, -7,  -2,  -9,  -10, -11, -6,  -1,  -8,  -3,  -4,
9959371c9d4SSatish Balay     8,   4,   0,   11,  1,   3,   2,   10,  6,   5,   7,   9,   -4,  -11, -9,  -1,  -8,  -12, -10, -5,  -3,  -7,  -2,  -6,  9,   11,  10,  6,   8,   7,   3,   5,   4,   0,   2,   1,   -3,  -1,  -2,  -6,  -4,  -5,  -9,  -7,  -8,  -12, -10, -11,
9969371c9d4SSatish Balay     10,  9,   11,  7,   6,   8,   4,   3,   5,   1,   0,   2,   -2,  -3,  -1,  -5,  -6,  -4,  -8,  -9,  -7,  -11, -12, -10, 11,  10,  9,   8,   7,   6,   5,   4,   3,   2,   1,   0,   -1,  -2,  -3,  -4,  -5,  -6,  -7,  -8,  -9,  -10, -11, -12,
9979371c9d4SSatish Balay     -12, -11, -10, -9,  -8,  -7,  -6,  -5,  -4,  -3,  -2,  -1,  0,   1,   2,   3,   4,   5,   6,   7,   8,   9,   10,  11,  -11, -10, -12, -8,  -7,  -9,  -5,  -4,  -6,  -2,  -1,  -3,  1,   2,   0,   4,   5,   3,   7,   8,   6,   10,  11,  9,
9989371c9d4SSatish Balay     -10, -12, -11, -7,  -9,  -8,  -4,  -6,  -5,  -1,  -3,  -2,  2,   0,   1,   5,   3,   4,   8,   6,   7,   11,  9,   10,  -9,  -5,  -1,  -12, -2,  -4,  -3,  -11, -7,  -6,  -8,  -10, 3,   10,  8,   0,   7,   11,  9,   4,   2,   6,   1,   5,
9999371c9d4SSatish Balay     -8,  -4,  -3,  -11, -1,  -6,  -2,  -10, -9,  -5,  -7,  -12, 4,   11,  6,   1,   8,   9,   10,  5,   0,   7,   2,   3,   -7,  -6,  -2,  -10, -3,  -5,  -1,  -12, -8,  -4,  -9,  -11, 5,   9,   7,   2,   6,   10,  11,  3,   1,   8,   0,   4,
10009371c9d4SSatish Balay     -3,  -8,  -4,  -6,  -11, -1,  -9,  -2,  -10, -12, -5,  -7,  6,   4,   11,  9,   1,   8,   0,   10,  5,   3,   7,   2,   -2,  -7,  -6,  -5,  -10, -3,  -8,  -1,  -12, -11, -4,  -9,  7,   5,   9,   10,  2,   6,   1,   11,  3,   4,   8,   0,
10019371c9d4SSatish Balay     -1,  -9,  -5,  -4,  -12, -2,  -7,  -3,  -11, -10, -6,  -8,  8,   3,   10,  11,  0,   7,   2,   9,   4,   5,   6,   1,   -6,  -2,  -7,  -3,  -5,  -10, -12, -8,  -1,  -9,  -11, -4,  9,   7,   5,   6,   10,  2,   3,   1,   11,  0,   4,   8,
10029371c9d4SSatish Balay     -5,  -1,  -9,  -2,  -4,  -12, -11, -7,  -3,  -8,  -10, -6,  10,  8,   3,   7,   11,  0,   4,   2,   9,   1,   5,   6,   -4,  -3,  -8,  -1,  -6,  -11, -10, -9,  -2,  -7,  -12, -5,  11,  6,   4,   8,   9,   1,   5,   0,   10,  2,   3,   7,
1003b5a892a1SMatthew G. Knepley   };
1004ef8b56bfSJed Brown   static const PetscInt hexMult[48 * 48] = {
1005b5a892a1SMatthew 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,
1006b5a892a1SMatthew 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,
1007b5a892a1SMatthew 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,
1008b5a892a1SMatthew 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,
1009b5a892a1SMatthew 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,
1010b5a892a1SMatthew 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,
1011b5a892a1SMatthew 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,
1012b5a892a1SMatthew 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,
1013b5a892a1SMatthew 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,
1014b5a892a1SMatthew 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,
1015b5a892a1SMatthew 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,
1016b5a892a1SMatthew 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,
1017b5a892a1SMatthew 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,
1018b5a892a1SMatthew 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,
1019b5a892a1SMatthew 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,
1020b5a892a1SMatthew 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,
1021b5a892a1SMatthew 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,
1022b5a892a1SMatthew 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,
1023b5a892a1SMatthew 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,
1024b5a892a1SMatthew 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,
1025b5a892a1SMatthew 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,
1026b5a892a1SMatthew 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,
1027b5a892a1SMatthew 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,
1028b5a892a1SMatthew 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,
1029b5a892a1SMatthew 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,
1030b5a892a1SMatthew 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,
1031b5a892a1SMatthew 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,
1032b5a892a1SMatthew 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,
1033b5a892a1SMatthew 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,
1034b5a892a1SMatthew 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,
1035b5a892a1SMatthew 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,
1036b5a892a1SMatthew 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,
1037b5a892a1SMatthew 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,
1038b5a892a1SMatthew 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,
1039b5a892a1SMatthew 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,
1040b5a892a1SMatthew 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,
1041b5a892a1SMatthew 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,
1042b5a892a1SMatthew 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,
1043b5a892a1SMatthew 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,
1044b5a892a1SMatthew 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,
1045b5a892a1SMatthew 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,
1046b5a892a1SMatthew 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,
1047b5a892a1SMatthew 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,
1048b5a892a1SMatthew 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,
1049b5a892a1SMatthew 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,
1050b5a892a1SMatthew 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,
1051b5a892a1SMatthew 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,
1052b5a892a1SMatthew 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,
1053b5a892a1SMatthew G. Knepley   };
1054ef8b56bfSJed Brown   static const PetscInt tripMult[12 * 12] = {
10559371c9d4SSatish Balay     1,  0,  2,  3,  5,  4,  -6, -4, -5, -2, -3, -1, 0,  2,  1,  4,  3,  5,  -5, -6, -4, -3, -1, -2, 2,  1,  0,  5,  4,  3,  -4, -5, -6, -1, -2, -3, 4,  3,  5,  0,  2,  1,  -3, -1, -2, -5, -6, -4,
10569371c9d4SSatish Balay     3,  5,  4,  1,  0,  2,  -2, -3, -1, -6, -4, -5, 5,  4,  3,  2,  1,  0,  -1, -2, -3, -4, -5, -6, -6, -5, -4, -3, -2, -1, 0,  1,  2,  3,  4,  5,  -4, -6, -5, -2, -1, -3, 1,  2,  0,  5,  3,  4,
10579371c9d4SSatish Balay     -5, -4, -6, -1, -3, -2, 2,  0,  1,  4,  5,  3,  -3, -2, -1, -6, -5, -4, 3,  4,  5,  0,  1,  2,  -1, -3, -2, -5, -4, -6, 4,  5,  3,  2,  0,  1,  -2, -1, -3, -4, -6, -5, 5,  3,  4,  1,  2,  0,
1058b5a892a1SMatthew G. Knepley   };
1059ef8b56bfSJed Brown   static const PetscInt ttriMult[12 * 12] = {
10609371c9d4SSatish Balay     0,  2,  1,  3,  5,  4,  -6, -4, -5, -3, -1, -2, 1,  0,  2,  4,  3,  5,  -5, -6, -4, -2, -3, -1, 2,  1,  0,  5,  4,  3,  -4, -5, -6, -1, -2, -3, 3,  5,  4,  0,  2,  1,  -3, -1, -2, -6, -4, -5,
10619371c9d4SSatish Balay     4,  3,  5,  1,  0,  2,  -2, -3, -1, -5, -6, -4, 5,  4,  3,  2,  1,  0,  -1, -2, -3, -4, -5, -6, -6, -5, -4, -3, -2, -1, 0,  1,  2,  3,  4,  5,  -5, -4, -6, -2, -1, -3, 1,  2,  0,  4,  5,  3,
10629371c9d4SSatish Balay     -4, -6, -5, -1, -3, -2, 2,  0,  1,  5,  3,  4,  -3, -2, -1, -6, -5, -4, 3,  4,  5,  0,  1,  2,  -2, -1, -3, -5, -4, -6, 4,  5,  3,  1,  2,  0,  -1, -3, -2, -4, -6, -5, 5,  3,  4,  2,  0,  1,
1063b5a892a1SMatthew G. Knepley   };
1064ef8b56bfSJed Brown   static const PetscInt tquadMult[16 * 16] = {
10659371c9d4SSatish Balay     0,  3,  2,  1,  4,  7,  6,  5,  -8, -5, -6, -7, -4, -1, -2, -3, 1,  0,  3,  2,  5,  4,  7,  6,  -7, -8, -5, -6, -3, -4, -1, -2, 2,  1,  0,  3,  6,  5,  4,  7,  -6, -7, -8, -5, -2, -3, -4, -1, 3, 2, 1, 0,
10669371c9d4SSatish Balay     7,  6,  5,  4,  -5, -6, -7, -8, -1, -2, -3, -4, 4,  7,  6,  5,  0,  3,  2,  1,  -4, -1, -2, -3, -8, -5, -6, -7, 5,  4,  7,  6,  1,  0,  3,  2,  -3, -4, -1, -2, -7, -8, -5, -6, 6,  5,  4,  7,  2, 1, 0, 3,
10679371c9d4SSatish Balay     -2, -3, -4, -1, -6, -7, -8, -5, 7,  6,  5,  4,  3,  2,  1,  0,  -1, -2, -3, -4, -5, -6, -7, -8, -8, -7, -6, -5, -4, -3, -2, -1, 0,  1,  2,  3,  4,  5,  6,  7,  -7, -6, -5, -8, -3, -2, -1, -4, 1, 2, 3, 0,
10689371c9d4SSatish Balay     5,  6,  7,  4,  -6, -5, -8, -7, -2, -1, -4, -3, 2,  3,  0,  1,  6,  7,  4,  5,  -5, -8, -7, -6, -1, -4, -3, -2, 3,  0,  1,  2,  7,  4,  5,  6,  -4, -3, -2, -1, -8, -7, -6, -5, 4,  5,  6,  7,  0, 1, 2, 3,
10699371c9d4SSatish Balay     -3, -2, -1, -4, -7, -6, -5, -8, 5,  6,  7,  4,  1,  2,  3,  0,  -2, -1, -4, -3, -6, -5, -8, -7, 6,  7,  4,  5,  2,  3,  0,  1,  -1, -4, -3, -2, -5, -8, -7, -6, 7,  4,  5,  6,  3,  0,  1,  2,
1070b5a892a1SMatthew G. Knepley   };
1071ef8b56bfSJed Brown   static const PetscInt pyrMult[8 * 8] = {
10729371c9d4SSatish Balay     0, 3, 2, 1, -4, -1, -2, -3, 1, 0, 3, 2, -3, -4, -1, -2, 2, 1, 0, 3, -2, -3, -4, -1, 3, 2, 1, 0, -1, -2, -3, -4, -4, -3, -2, -1, 0, 1, 2, 3, -3, -2, -1, -4, 1, 2, 3, 0, -2, -1, -4, -3, 2, 3, 0, 1, -1, -4, -3, -2, 3, 0, 1, 2,
1073b5a892a1SMatthew G. Knepley   };
1074b5a892a1SMatthew G. Knepley   switch (ct) {
1075d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
1076d71ae5a4SJacob Faibussowitsch     return 0;
1077b5a892a1SMatthew G. Knepley   case DM_POLYTOPE_SEGMENT:
1078d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
1079d71ae5a4SJacob Faibussowitsch     return segMult[(o1 + 1) * 2 + o2 + 1];
1080d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
1081d71ae5a4SJacob Faibussowitsch     return triMult[(o1 + 3) * 6 + o2 + 3];
1082d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
1083d71ae5a4SJacob Faibussowitsch     return quadMult[(o1 + 4) * 8 + o2 + 4];
1084d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
1085d71ae5a4SJacob Faibussowitsch     return tsegMult[(o1 + 2) * 4 + o2 + 2];
1086d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
1087d71ae5a4SJacob Faibussowitsch     return tetMult[(o1 + 12) * 24 + o2 + 12];
1088d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
1089d71ae5a4SJacob Faibussowitsch     return hexMult[(o1 + 24) * 48 + o2 + 24];
1090d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
1091d71ae5a4SJacob Faibussowitsch     return tripMult[(o1 + 6) * 12 + o2 + 6];
1092d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
1093d71ae5a4SJacob Faibussowitsch     return ttriMult[(o1 + 6) * 12 + o2 + 6];
1094d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
1095d71ae5a4SJacob Faibussowitsch     return tquadMult[(o1 + 8) * 16 + o2 + 8];
1096d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
1097d71ae5a4SJacob Faibussowitsch     return pyrMult[(o1 + 4) * 8 + o2 + 4];
1098d71ae5a4SJacob Faibussowitsch   default:
1099d71ae5a4SJacob Faibussowitsch     return 0;
1100b5a892a1SMatthew G. Knepley   }
1101b5a892a1SMatthew G. Knepley }
1102b5a892a1SMatthew G. Knepley 
1103b5a892a1SMatthew G. Knepley /* This is orientation o1 acting on orientation o2^{-1} */
DMPolytopeTypeComposeOrientationInv(DMPolytopeType ct,PetscInt o1,PetscInt o2)1104d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeComposeOrientationInv(DMPolytopeType ct, PetscInt o1, PetscInt o2)
1105d71ae5a4SJacob Faibussowitsch {
1106ef8b56bfSJed Brown   static const PetscInt triInv[6]    = {-3, -2, -1, 0, 2, 1};
1107ef8b56bfSJed Brown   static const PetscInt quadInv[8]   = {-4, -3, -2, -1, 0, 3, 2, 1};
1108ef8b56bfSJed 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};
11099371c9d4SSatish Balay   static const PetscInt hexInv[48]   = {-17, -18, -20, -19, -22, -21, -23, -24, -15, -16, -14, -13, -11, -12, -10, -9, -8, -5, -6, -7, -4, -3, -2, -1, 0, 3, 2, 1, 6, 5, 4, 9, 8, 7, 10, 11, 12, 13, 14, 15, 17, 16, 19, 18, 21, 20, 23, 22};
1110ef8b56bfSJed Brown   static const PetscInt tripInv[12]  = {-5, -6, -4, -3, -2, -1, 0, 2, 1, 3, 4, 5};
1111ef8b56bfSJed Brown   static const PetscInt ttriInv[12]  = {-6, -5, -4, -3, -2, -1, 0, 2, 1, 3, 5, 4};
1112ef8b56bfSJed Brown   static const PetscInt tquadInv[16] = {-8, -7, -6, -5, -4, -3, -2, -1, 0, 3, 2, 1, 4, 7, 6, 5};
1113ef8b56bfSJed Brown   static const PetscInt pyrInv[8]    = {-4, -3, -2, -1, 0, 3, 2, 1};
1114b5a892a1SMatthew G. Knepley   switch (ct) {
1115d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT:
1116d71ae5a4SJacob Faibussowitsch     return 0;
1117b5a892a1SMatthew G. Knepley   case DM_POLYTOPE_SEGMENT:
1118d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_POINT_PRISM_TENSOR:
1119d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, o2);
1120d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRIANGLE:
1121d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, triInv[o2 + 3]);
1122d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUADRILATERAL:
1123d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, quadInv[o2 + 4]);
1124d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_SEG_PRISM_TENSOR:
1125d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, o2);
1126d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TETRAHEDRON:
1127d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, tetInv[o2 + 12]);
1128d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_HEXAHEDRON:
1129d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, hexInv[o2 + 24]);
1130d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM:
1131d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, tripInv[o2 + 6]);
1132d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_TRI_PRISM_TENSOR:
1133d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, ttriInv[o2 + 6]);
1134d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_QUAD_PRISM_TENSOR:
1135d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, tquadInv[o2 + 8]);
1136d71ae5a4SJacob Faibussowitsch   case DM_POLYTOPE_PYRAMID:
1137d71ae5a4SJacob Faibussowitsch     return DMPolytopeTypeComposeOrientation(ct, o1, pyrInv[o2 + 4]);
1138d71ae5a4SJacob Faibussowitsch   default:
1139d71ae5a4SJacob Faibussowitsch     return 0;
1140b5a892a1SMatthew G. Knepley   }
1141b5a892a1SMatthew G. Knepley }
1142b5a892a1SMatthew G. Knepley 
1143b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeMatchOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *);
1144b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeMatchVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *);
1145b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeGetOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *);
1146b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeGetVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *);
1147012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeInCellTest(DMPolytopeType, const PetscReal[], PetscBool *);
1148