1e1589f56SBarry Smith /* 2e1589f56SBarry Smith Objects to manage the interactions between the mesh data structures and the algebraic objects 3e1589f56SBarry Smith */ 4a4963045SJacob Faibussowitsch #pragma once 52c8e378dSBarry Smith #include <petscmat.h> 61e25c274SJed Brown #include <petscdmtypes.h> 7decb47aaSMatthew G. Knepley #include <petscfetypes.h> 82764a2aaSMatthew G. Knepley #include <petscdstypes.h> 9c58f1c22SToby Isaac #include <petscdmlabel.h> 1007218a29SMatthew G. Knepley #include <petscdt.h> 11e1589f56SBarry Smith 12ac09b921SBarry Smith /* SUBMANSEC = DM */ 13ac09b921SBarry Smith 14607a6623SBarry Smith PETSC_EXTERN PetscErrorCode DMInitializePackage(void); 15e1589f56SBarry Smith 16014dd563SJed Brown PETSC_EXTERN PetscClassId DM_CLASSID; 17e1589f56SBarry Smith 18528c63e0SDave May #define DMLOCATEPOINT_POINT_NOT_FOUND -367 19528c63e0SDave May 2076bdecfbSBarry Smith /*J 2187497f52SBarry Smith DMType - String with the name of a PETSc `DM` 22e1589f56SBarry Smith 23e1589f56SBarry Smith Level: beginner 24e1589f56SBarry Smith 251cc06b55SBarry Smith .seealso: [](ch_dmbase), `DMSetType()`, `DMCreate()`, `DM` 2676bdecfbSBarry Smith J*/ 2719fd82e9SBarry Smith typedef const char *DMType; 28e1589f56SBarry Smith #define DMDA "da" 29e1589f56SBarry Smith #define DMCOMPOSITE "composite" 30e1589f56SBarry Smith #define DMSLICED "sliced" 31fe1899a2SJed Brown #define DMSHELL "shell" 32ab7f58a0SBarry Smith #define DMPLEX "plex" 338ac4e037SJed Brown #define DMREDUNDANT "redundant" 343a19ef87SMatthew G Knepley #define DMPATCH "patch" 351d72bce8STim Tautges #define DMMOAB "moab" 365f2c45f1SShri Abhyankar #define DMNETWORK "network" 37ef51cf95SToby Isaac #define DMFOREST "forest" 38db4d5e8cSToby Isaac #define DMP4EST "p4est" 39db4d5e8cSToby Isaac #define DMP8EST "p8est" 402fd35b1fSDave May #define DMSWARM "swarm" 41d852a638SPatrick Sanan #define DMPRODUCT "product" 42a3101111SPatrick Sanan #define DMSTAG "stag" 43e1589f56SBarry Smith 44bff4a2f0SMatthew G. Knepley PETSC_EXTERN const char *const DMBoundaryTypes[]; 4540967b3bSMatthew G. Knepley PETSC_EXTERN const char *const DMBoundaryConditionTypes[]; 46863027abSJed Brown PETSC_EXTERN const char *const DMBlockingTypes[]; 47140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList DMList; 48c0517cd5SMatthew G. Knepley PETSC_EXTERN DMGeneratorFunctionList DMGenerateList; 497625649eSMatthew G. Knepley PETSC_EXTERN PetscFunctionList DMGeomModelList; 50014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreate(MPI_Comm, DM *); 5138221697SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMClone(DM, DM *); 5219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetType(DM, DMType); 5319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMGetType(DM, DMType *); 54bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode DMRegister(const char[], PetscErrorCode (*)(DM)); 55014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRegisterDestroy(void); 56e1589f56SBarry Smith 57014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMView(DM, PetscViewer); 58014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLoad(DM, PetscViewer); 59014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMDestroy(DM *); 60014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM, Vec *); 61014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM, Vec *); 62014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM, Vec *); 63014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM, Vec *); 64014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM, Vec *); 65014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM, Vec *); 66014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM); 6750eeb1caSToby Isaac PETSC_EXTERN PetscErrorCode DMClearLocalVectors(DM); 68974ca4ecSStefano Zampini PETSC_EXTERN PetscErrorCode DMClearNamedGlobalVectors(DM); 69974ca4ecSStefano Zampini PETSC_EXTERN PetscErrorCode DMClearNamedLocalVectors(DM); 70e77ac854SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasNamedGlobalVector(DM, const char *, PetscBool *); 71014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM, const char *, Vec *); 72014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM, const char *, Vec *); 73e77ac854SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasNamedLocalVector(DM, const char *, PetscBool *); 742348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM, const char *, Vec *); 752348bcf4SPeter Brune PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM, const char *, Vec *); 76014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM, ISLocalToGlobalMapping *); 77014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM, PetscInt *, char ***, IS **); 78014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM, PetscInt *); 79b412c318SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateColoring(DM, ISColoringType, ISColoring *); 80b412c318SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM, Mat *); 81aa0f6e3cSJed Brown PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateSkip(DM, PetscBool); 82014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM, PetscBool); 83b06ff27eSHong Zhang PETSC_EXTERN PetscErrorCode DMSetMatrixStructureOnly(DM, PetscBool); 84863027abSJed Brown PETSC_EXTERN PetscErrorCode DMSetBlockingType(DM, DMBlockingType); 85863027abSJed Brown PETSC_EXTERN PetscErrorCode DMGetBlockingType(DM, DMBlockingType *); 86014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM, DM, Mat *, Vec *); 873ad4599aSBarry Smith PETSC_EXTERN PetscErrorCode DMCreateRestriction(DM, DM, Mat *); 886dbf9973SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMCreateInjection(DM, DM, Mat *); 89bd041c0cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateMassMatrix(DM, DM, Mat *); 90*8e9849d2SStefano Zampini PETSC_EXTERN PetscErrorCode DMCreateMassMatrixLumped(DM, Vec *, Vec *); 9169291d52SBarry Smith PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM, PetscInt, MPI_Datatype, void *); 9269291d52SBarry Smith PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM, PetscInt, MPI_Datatype, void *); 93014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefine(DM, MPI_Comm, DM *); 94014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsen(DM, MPI_Comm, DM *); 95a8fb8f29SToby Isaac PETSC_EXTERN PetscErrorCode DMGetCoarseDM(DM, DM *); 96a8fb8f29SToby Isaac PETSC_EXTERN PetscErrorCode DMSetCoarseDM(DM, DM); 9788bdff64SToby Isaac PETSC_EXTERN PetscErrorCode DMGetFineDM(DM, DM *); 9888bdff64SToby Isaac PETSC_EXTERN PetscErrorCode DMSetFineDM(DM, DM); 99014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM, PetscInt, DM[]); 100014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM, PetscInt, DM[]); 101014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, Vec, Mat, DM, void *), void *); 102dc822a44SJed Brown PETSC_EXTERN PetscErrorCode DMCoarsenHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, Vec, Mat, DM, void *), void *); 103014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, DM, void *), void *); 1043d8e3701SJed Brown PETSC_EXTERN PetscErrorCode DMRefineHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, DM, void *), void *); 105014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMRestrict(DM, Mat, Vec, Mat, DM); 106014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMInterpolate(DM, Mat, DM); 1071f3379b2SToby Isaac PETSC_EXTERN PetscErrorCode DMInterpolateSolution(DM, DM, Mat, Vec, Vec); 108d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMExtrude(DM, PetscInt, DM *); 109014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM); 110fe2efc57SMark PETSC_EXTERN PetscErrorCode DMViewFromOptions(DM, PetscObject, const char[]); 111ca266f36SBarry Smith 112c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerate(DM, const char[], PetscBool, DM *); 1139fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMGenerateRegister(const char[], PetscErrorCode (*)(DM, PetscBool, DM *), PetscErrorCode (*)(DM, PetscReal *, DM *), PetscErrorCode (*)(DM, Vec, DMLabel, DMLabel, DM *), PetscInt); 114c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerateRegisterAll(void); 115c0517cd5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGenerateRegisterDestroy(void); 1167625649eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGeomModelRegister(const char[], PetscErrorCode (*)(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[])); 1177625649eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGeomModelRegisterAll(void); 1187625649eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGeomModelRegisterDestroy(void); 119a1b0c543SToby Isaac PETSC_EXTERN PetscErrorCode DMAdaptLabel(DM, DMLabel, DM *); 1209fe9e680SJoe Wallwork PETSC_EXTERN PetscErrorCode DMAdaptMetric(DM, Vec, DMLabel, DMLabel, DM *); 121df0b854cSToby Isaac 122014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetUp(DM); 123014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM, DM, Mat, Vec *); 124edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMDACreateAggregates()", ) PetscErrorCode DMCreateAggregates(DM, DM, Mat *); 125baf369e7SPeter Brune PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM, PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), void *); 126d4d07f1eSToby Isaac PETSC_EXTERN PetscErrorCode DMLocalToGlobalHookAdd(DM, PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), void *); 12701729b5cSPatrick Sanan PETSC_EXTERN PetscErrorCode DMGlobalToLocal(DM, Vec, InsertMode, Vec); 128014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM, Vec, InsertMode, Vec); 129014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM, Vec, InsertMode, Vec); 13001729b5cSPatrick Sanan PETSC_EXTERN PetscErrorCode DMLocalToGlobal(DM, Vec, InsertMode, Vec); 131014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM, Vec, InsertMode, Vec); 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM, Vec, InsertMode, Vec); 133d78e899eSRichard Tran Mills PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM, Vec, InsertMode, Vec); 134d78e899eSRichard Tran Mills PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM, Vec, InsertMode, Vec); 13519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMConvert(DM, DMType, DM *); 136e1589f56SBarry Smith 137c73cfb54SMatthew G. Knepley /* Topology support */ 138c73cfb54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDimension(DM, PetscInt *); 139c73cfb54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetDimension(DM, PetscInt); 140793f3fe5SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDimPoints(DM, PetscInt, PetscInt *, PetscInt *); 1418e4ac7eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetUseNatural(DM, PetscBool *); 1428e4ac7eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetUseNatural(DM, PetscBool); 1436858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNeighbors(DM, PetscInt *, const PetscMPIInt **); 14446e270d4SMatthew G. Knepley 14546e270d4SMatthew G. Knepley /* Coordinate support */ 1466636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM, DM *); 1471cfe2091SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM, DM); 1486858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinateDM(DM, DM *); 1496858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinateDM(DM, DM); 15046e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateDim(DM, PetscInt *); 15146e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateDim(DM, PetscInt); 152e8abe2deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM, PetscSection *); 15346e270d4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM, PetscInt, PetscSection); 1546858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinateSection(DM, PetscSection *); 1556858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinateSection(DM, PetscInt, PetscSection); 1566636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM, Vec *); 1576636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM, Vec); 1586858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinates(DM, Vec *); 1596858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinates(DM, Vec); 16081e9a530SVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalSetUp(DM); 1616858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM, Vec *); 16281e9a530SVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalNoncollective(DM, Vec *); 1632db98f8dSVaclav Hapla PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalTuple(DM, IS, PetscSection *, Vec *); 1646636e97aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM, Vec); 1656858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM); 1666858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocal(DM, Vec *); 1676858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM, Vec *); 1686858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCellCoordinatesLocal(DM, Vec); 1696858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCoordinateField(DM, DMField *); 1706858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateField(DM, DMField); 1716858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetLocalBoundingBox(DM, PetscReal[], PetscReal[]); 1726858538eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetBoundingBox(DM, PetscReal[], PetscReal[]); 173e44f6aebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoordinateDisc(DM, PetscFE, PetscBool); 17462a38674SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocatePoints(DM, Vec, DMPointLocationType, PetscSF *); 1757625649eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSnapToGeomModel(DM, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]); 1767625649eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetSnapToGeomModel(DM, const char[]); 1776858538eSMatthew G. Knepley 1786858538eSMatthew G. Knepley /* Periodicity support */ 1794fb89dddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM, const PetscReal *[], const PetscReal *[], const PetscReal *[]); 1804fb89dddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM, const PetscReal[], const PetscReal[], const PetscReal[]); 181e907e85cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate(DM, const PetscScalar[], PetscBool, PetscScalar[]); 1822e17dfb7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMLocalizeCoordinates(DM); 18336447a5eSToby Isaac PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalized(DM, PetscBool *); 1848f700142SStefano Zampini PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalizedLocal(DM, PetscBool *); 185a4b8866cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetSparseLocalize(DM, PetscBool *); 186a4b8866cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetSparseLocalize(DM, PetscBool); 1876636e97aSMatthew G Knepley 1885dbd56e3SPeter Brune /* block hook interface */ 189be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, VecScatter, VecScatter, DM, void *), void *); 190b3a6b972SJed Brown PETSC_EXTERN PetscErrorCode DMSubDomainHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, VecScatter, VecScatter, DM, void *), void *); 191be081cd6SPeter Brune PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM, VecScatter, VecScatter, DM); 1925dbd56e3SPeter Brune 193014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM, const char[]); 19431697293SDave May PETSC_EXTERN PetscErrorCode DMAppendOptionsPrefix(DM, const char[]); 19531697293SDave May PETSC_EXTERN PetscErrorCode DMGetOptionsPrefix(DM, const char *[]); 19619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetVecType(DM, VecType); 197c0dedaeaSBarry Smith PETSC_EXTERN PetscErrorCode DMGetVecType(DM, VecType *); 19819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode DMSetMatType(DM, MatType); 199c0dedaeaSBarry Smith PETSC_EXTERN PetscErrorCode DMGetMatType(DM, MatType *); 2008f1509bcSBarry Smith PETSC_EXTERN PetscErrorCode DMSetISColoringType(DM, ISColoringType); 2018f1509bcSBarry Smith PETSC_EXTERN PetscErrorCode DMGetISColoringType(DM, ISColoringType *); 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM, void *); 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM, PetscErrorCode (*)(void **)); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM, void *); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM, PetscErrorCode (*)(DM, Vec, Vec)); 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM, PetscBool *); 207b0ae01b7SPeter Brune PETSC_EXTERN PetscErrorCode DMHasColoring(DM, PetscBool *); 2083ad4599aSBarry Smith PETSC_EXTERN PetscErrorCode DMHasCreateRestriction(DM, PetscBool *); 209a7058e45SLawrence Mitchell PETSC_EXTERN PetscErrorCode DMHasCreateInjection(DM, PetscBool *); 210014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM, Vec, Vec); 21193d92d96SBarry Smith 21237bc7515SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, const PetscInt[], IS *, DM *); 2132adcc780SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSuperDM(DM[], PetscInt, IS **, DM *); 214dd072f5fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSectionSubDM(DM, PetscInt, const PetscInt[], const PetscInt[], const PetscInt[], IS *, DM *); 215792b654fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSectionSuperDM(DM[], PetscInt, IS **, DM *); 21616621825SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM, PetscInt *, char ***, IS **, DM **); 2178d4ac253SDmitry Karpeev PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM, PetscInt *, char ***, IS **, IS **, DM **); 218e30e807fSPeter Brune PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM, PetscInt, DM *, VecScatter **, VecScatter **, VecScatter **); 219e7c4fc90SDmitry Karpeev 220014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM, PetscInt *); 221fef3a512SBarry Smith PETSC_EXTERN PetscErrorCode DMSetRefineLevel(DM, PetscInt); 222014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM, PetscInt *); 2239a64c4a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetCoarsenLevel(DM, PetscInt); 224014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMFinalizePackage(void); 225e1589f56SBarry Smith 2265f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM *); 2275f1ad066SMatthew G Knepley PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM); 228c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM *); 229c688c046SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM); 230531c7667SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringUseDM(Mat, MatFDColoring); 2315f1ad066SMatthew G Knepley 232e1589f56SBarry Smith typedef struct NLF_DAAD *NLF; 233e1589f56SBarry Smith 234bc2bf880SBarry Smith #define DM_FILE_CLASSID 1211221 2357da65231SMatthew G Knepley 2367da65231SMatthew G Knepley /* FEM support */ 2372ecf6ec3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPrintCellIndices(PetscInt, const char[], PetscInt, const PetscInt[]); 238014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char[], PetscInt, const PetscScalar[]); 2395962854dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPrintCellVectorReal(PetscInt, const char[], PetscInt, const PetscReal[]); 240014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char[], PetscInt, PetscInt, const PetscScalar[]); 2416113b454SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char[], PetscReal, Vec); 2427da65231SMatthew G Knepley 2438cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *)); 2448cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *)); 2458cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *)); 2468cda7954SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *)); 247f9d4088aSMatthew G. Knepley 248061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMGetSection(DM, PetscSection *); /* Use DMGetLocalSection() in new code (since v3.12) */ 249061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMSetSection(DM, PetscSection); /* Use DMSetLocalSection() in new code (since v3.12) */ 250061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMGetLocalSection(DM, PetscSection *); 251061576a5SJed Brown PETSC_EXTERN PetscErrorCode DMSetLocalSection(DM, PetscSection); 252e87a4003SBarry Smith PETSC_EXTERN PetscErrorCode DMGetGlobalSection(DM, PetscSection *); 253e87a4003SBarry Smith PETSC_EXTERN PetscErrorCode DMSetGlobalSection(DM, PetscSection); 254adc21957SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateSectionPermutation(DM, IS *, PetscBT *); 255adc21957SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMReorderSectionGetDefault(DM, DMReorderDefaultFlag *); 256adc21957SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMReorderSectionSetDefault(DM, DMReorderDefaultFlag); 257adc21957SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMReorderSectionGetType(DM, MatOrderingType *); 258adc21957SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMReorderSectionSetType(DM, MatOrderingType); 259d2b2dc1eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMUseTensorOrder(DM, PetscBool); 260edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMGetSection()", ) PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *s) 261d71ae5a4SJacob Faibussowitsch { 2629371c9d4SSatish Balay return DMGetSection(dm, s); 2639371c9d4SSatish Balay } 264edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMSetSection()", ) PetscErrorCode DMSetDefaultSection(DM dm, PetscSection s) 265d71ae5a4SJacob Faibussowitsch { 2669371c9d4SSatish Balay return DMSetSection(dm, s); 2679371c9d4SSatish Balay } 268edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMGetGlobalSection()", ) PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *s) 269d71ae5a4SJacob Faibussowitsch { 2709371c9d4SSatish Balay return DMGetGlobalSection(dm, s); 2719371c9d4SSatish Balay } 272edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 9, 0, "DMSetGlobalSection()", ) PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection s) 273d71ae5a4SJacob Faibussowitsch { 2749371c9d4SSatish Balay return DMSetGlobalSection(dm, s); 2759371c9d4SSatish Balay } 276e87a4003SBarry Smith 2771bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMGetSectionSF(DM, PetscSF *); 2781bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMSetSectionSF(DM, PetscSF); 2791bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMCreateSectionSF(DM, PetscSection, PetscSection); 280edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMGetSectionSF()", ) PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *s) 281d71ae5a4SJacob Faibussowitsch { 2829371c9d4SSatish Balay return DMGetSectionSF(dm, s); 2839371c9d4SSatish Balay } 284edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMSetSectionSF()", ) PetscErrorCode DMSetDefaultSF(DM dm, PetscSF s) 285d71ae5a4SJacob Faibussowitsch { 2869371c9d4SSatish Balay return DMSetSectionSF(dm, s); 2879371c9d4SSatish Balay } 288edd03b47SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION(3, 12, 0, "DMCreateSectionSF()", ) PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection l, PetscSection g) 289d71ae5a4SJacob Faibussowitsch { 2909371c9d4SSatish Balay return DMCreateSectionSF(dm, l, g); 2919371c9d4SSatish Balay } 2921bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *); 2931bb6d2a8SBarry Smith PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF); 2944f37162bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNaturalSF(DM, PetscSF *); 2954f37162bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetNaturalSF(DM, PetscSF); 2961bb6d2a8SBarry Smith 29779769bd5SJed Brown PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *, Vec *); 29879769bd5SJed Brown PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat, Vec); 299e87a4003SBarry Smith 30014f150ffSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *); 301cdb7a50dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *); 302cdb7a50dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal); 303b2033f5dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char[], PetscInt, PetscReal *); 304b2033f5dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetOutputSequenceLength(DM, PetscViewer, const char[], PetscInt *); 30514f150ffSMatthew G. Knepley 306af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *); 307af122d2aSMatthew G Knepley PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt); 30844a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, DMLabel *, PetscObject *); 30944a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, DMLabel, PetscObject); 31044a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAddField(DM, DMLabel, PetscObject); 311e0b68406SMatthew Knepley PETSC_EXTERN PetscErrorCode DMSetFieldAvoidTensor(DM, PetscInt, PetscBool); 312e0b68406SMatthew Knepley PETSC_EXTERN PetscErrorCode DMGetFieldAvoidTensor(DM, PetscInt, PetscBool *); 31344a7f3ddSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMClearFields(DM); 314e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyFields(DM, DM); 31534aa8a36SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAdjacency(DM, PetscInt, PetscBool *, PetscBool *); 31634aa8a36SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetAdjacency(DM, PetscInt, PetscBool, PetscBool); 317b0441da4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetBasicAdjacency(DM, PetscBool *, PetscBool *); 318b0441da4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetBasicAdjacency(DM, PetscBool, PetscBool); 319e5e52638SMatthew G. Knepley 320e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNumDS(DM, PetscInt *); 321e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *); 32207218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetCellDS(DM, PetscInt, PetscDS *, PetscDS *); 32307218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetRegionDS(DM, DMLabel, IS *, PetscDS *, PetscDS *); 32407218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetRegionDS(DM, DMLabel, IS, PetscDS, PetscDS); 32507218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetRegionNumDS(DM, PetscInt, DMLabel *, IS *, PetscDS *, PetscDS *); 32607218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetRegionNumDS(DM, PetscInt, DMLabel, IS, PetscDS, PetscDS); 3271d3af9e0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMFindRegionNum(DM, PetscDS, PetscInt *); 3282df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateFEDefault(DM, PetscInt, const char[], PetscInt, PetscFE *); 329e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCreateDS(DM); 330e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMClearDS(DM); 331e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyDS(DM, DM); 332e5e52638SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyDisc(DM, DM); 333f2cacb80SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMComputeExactSolution(DM, PetscReal, Vec, Vec); 3349a2a23afSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetNumAuxiliaryVec(DM, PetscInt *); 335ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAuxiliaryVec(DM, DMLabel, PetscInt, PetscInt, Vec *); 336ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetAuxiliaryVec(DM, DMLabel, PetscInt, PetscInt, Vec); 337ac17215fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMGetAuxiliaryLabels(DM, DMLabel[], PetscInt[], PetscInt[]); 3389a2a23afSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyAuxiliaryVec(DM, DM); 339e4d5475eSStefano Zampini PETSC_EXTERN PetscErrorCode DMClearAuxiliaryVec(DM); 340af122d2aSMatthew G Knepley 3414267b1a3SMatthew G. Knepley /*MC 3424267b1a3SMatthew G. Knepley DMInterpolationInfo - Structure for holding information about interpolation on a mesh 3434267b1a3SMatthew G. Knepley 3444267b1a3SMatthew G. Knepley Synopsis: 3454267b1a3SMatthew G. Knepley comm - The communicator 3464267b1a3SMatthew G. Knepley dim - The spatial dimension of points 3474267b1a3SMatthew G. Knepley nInput - The number of input points 3484267b1a3SMatthew G. Knepley points - The input point coordinates 3494267b1a3SMatthew G. Knepley cells - The cell containing each point 3504267b1a3SMatthew G. Knepley n - The number of local points 3514267b1a3SMatthew G. Knepley coords - The point coordinates 3524267b1a3SMatthew G. Knepley dof - The number of components to interpolate 3534267b1a3SMatthew G. Knepley 35416a05f60SBarry Smith Level: intermediate 35516a05f60SBarry Smith 356af27ebaaSBarry Smith .seealso: [](ch_dmbase), `DM`, `DMInterpolationCreate()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()` 3574267b1a3SMatthew G. Knepley M*/ 358e87bb0d3SMatthew G Knepley struct _DMInterpolationInfo { 359e87bb0d3SMatthew G Knepley MPI_Comm comm; 360e87bb0d3SMatthew G Knepley PetscInt dim; /*1 The spatial dimension of points */ 361e87bb0d3SMatthew G Knepley PetscInt nInput; /* The number of input points */ 362e87bb0d3SMatthew G Knepley PetscReal *points; /* The input point coordinates */ 363e87bb0d3SMatthew G Knepley PetscInt *cells; /* The cell containing each point */ 364e87bb0d3SMatthew G Knepley PetscInt n; /* The number of local points */ 365e87bb0d3SMatthew G Knepley Vec coords; /* The point coordinates */ 366e87bb0d3SMatthew G Knepley PetscInt dof; /* The number of components to interpolate */ 367e87bb0d3SMatthew G Knepley }; 368e87bb0d3SMatthew G Knepley typedef struct _DMInterpolationInfo *DMInterpolationInfo; 369e87bb0d3SMatthew G Knepley 37094b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *); 37194b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt); 37294b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *); 37394b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt); 37494b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *); 37594b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]); 37652aa1562SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool, PetscBool); 37794b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *); 37894b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *); 37994b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *); 38094b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec); 38194b4b8a8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *); 382c58f1c22SToby Isaac 383c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMCreateLabel(DM, const char[]); 3844bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreateLabelAtIndex(DM, PetscInt, const char[]); 385c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelValue(DM, const char[], PetscInt, PetscInt *); 386c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMSetLabelValue(DM, const char[], PetscInt, PetscInt); 387c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMClearLabelValue(DM, const char[], PetscInt, PetscInt); 388c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelSize(DM, const char[], PetscInt *); 389c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelIdIS(DM, const char[], IS *); 390c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetStratumSize(DM, const char[], PetscInt, PetscInt *); 391c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetStratumIS(DM, const char[], PetscInt, IS *); 3924de306b1SToby Isaac PETSC_EXTERN PetscErrorCode DMSetStratumIS(DM, const char[], PetscInt, IS); 393c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMClearLabelStratum(DM, const char[], PetscInt); 394c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelOutput(DM, const char[], PetscBool *); 395c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMSetLabelOutput(DM, const char[], PetscBool); 3965f15299fSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMGetFirstLabeledPoint(DM, DM, DMLabel, PetscInt, const PetscInt *, PetscInt, PetscInt *, PetscDS *); 397c58f1c22SToby Isaac 3982cbb9b06SVaclav Hapla /*E 399af27ebaaSBarry Smith DMCopyLabelsMode - Determines how `DMCopyLabels()` behaves when there is a `DMLabel` in the source and destination `DM`s with the same name 4002cbb9b06SVaclav Hapla 40116a05f60SBarry Smith Values: 40216a05f60SBarry Smith + `DM_COPY_LABELS_REPLACE` - replace label in destination by label from source 40316a05f60SBarry Smith . `DM_COPY_LABELS_KEEP` - keep destination label 404af27ebaaSBarry Smith - `DM_COPY_LABELS_FAIL` - generate an error 40516a05f60SBarry Smith 4062cbb9b06SVaclav Hapla Level: advanced 4072cbb9b06SVaclav Hapla 408af27ebaaSBarry Smith .seealso: [](ch_dmbase), `DMLabel`, `DM`, `DMCompareLabels()`, `DMRemoveLabel()` 4092cbb9b06SVaclav Hapla E*/ 4109371c9d4SSatish Balay typedef enum { 4119371c9d4SSatish Balay DM_COPY_LABELS_REPLACE, 4129371c9d4SSatish Balay DM_COPY_LABELS_KEEP, 4139371c9d4SSatish Balay DM_COPY_LABELS_FAIL 4149371c9d4SSatish Balay } DMCopyLabelsMode; 4152cbb9b06SVaclav Hapla PETSC_EXTERN const char *const DMCopyLabelsModes[]; 4162cbb9b06SVaclav Hapla 417c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetNumLabels(DM, PetscInt *); 418c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelName(DM, PetscInt, const char **); 419c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMHasLabel(DM, const char[], PetscBool *); 420c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabel(DM, const char *, DMLabel *); 4214a7ee7d0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSetLabel(DM, DMLabel); 422c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMGetLabelByNum(DM, PetscInt, DMLabel *); 423c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMAddLabel(DM, DMLabel); 424c58f1c22SToby Isaac PETSC_EXTERN PetscErrorCode DMRemoveLabel(DM, const char[], DMLabel *); 425306894acSVaclav Hapla PETSC_EXTERN PetscErrorCode DMRemoveLabelBySelf(DM, DMLabel *, PetscBool); 4262cbb9b06SVaclav Hapla PETSC_EXTERN PetscErrorCode DMCopyLabels(DM, DM, PetscCopyMode, PetscBool, DMCopyLabelsMode emode); 427609dae6eSVaclav Hapla PETSC_EXTERN PetscErrorCode DMCompareLabels(DM, DM, PetscBool *, char **); 428c58f1c22SToby Isaac 42945480ffeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAddBoundary(DM, DMBoundaryConditionType, const char[], DMLabel, PetscInt, const PetscInt[], PetscInt, PetscInt, const PetscInt[], void (*)(void), void (*)(void), void *, PetscInt *); 430a6ba4734SToby Isaac PETSC_EXTERN PetscErrorCode DMIsBoundaryPoint(DM, PetscInt, PetscBool *); 4314d6f44ffSToby Isaac 4320709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMProjectFunction(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec); 4330709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMProjectFunctionLocal(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec); 4342c53366bSMatthew 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); 4351c531cf8SMatthew 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); 436191494d9SMatthew 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); 437d29d7c6eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectFieldLabel(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Vec, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode, Vec); 4381c531cf8SMatthew 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); 439ece3a9fcSMatthew 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); 4400709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2Diff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *); 441b698f381SToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2GradientDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal[], PetscReal *); 4421189c1efSToby Isaac PETSC_EXTERN PetscErrorCode DMComputeL2FieldDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *); 4432e4af2aeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMComputeError(DM, Vec, PetscReal[], Vec *); 444ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMHasBasisTransform(DM, PetscBool *); 445ca3d3a14SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCopyTransform(DM, DM); 4468320bc6fSPatrick Sanan 4478320bc6fSPatrick Sanan PETSC_EXTERN PetscErrorCode DMGetCompatibility(DM, DM, PetscBool *, PetscBool *); 448c0f0dcc3SMatthew G. Knepley 449c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitorSet(DM, PetscErrorCode (*)(DM, void *), void *, PetscErrorCode (*)(void **)); 450c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitorCancel(DM); 451c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitorSetFromOptions(DM, const char[], const char[], const char[], PetscErrorCode (*)(DM, void *), PetscErrorCode (*)(DM, PetscViewerAndFormat *), PetscBool *); 452c0f0dcc3SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMMonitor(DM); 453c0f0dcc3SMatthew G. Knepley 454c306944fSJed Brown static inline PetscBool DMPolytopeTypeIsHybrid(DMPolytopeType ct) 455c306944fSJed Brown { 456c306944fSJed Brown switch (ct) { 457c306944fSJed Brown case DM_POLYTOPE_POINT_PRISM_TENSOR: 458c306944fSJed Brown case DM_POLYTOPE_SEG_PRISM_TENSOR: 459c306944fSJed Brown case DM_POLYTOPE_TRI_PRISM_TENSOR: 460c306944fSJed Brown case DM_POLYTOPE_QUAD_PRISM_TENSOR: 461c306944fSJed Brown return PETSC_TRUE; 462c306944fSJed Brown default: 463c306944fSJed Brown return PETSC_FALSE; 464c306944fSJed Brown } 465c306944fSJed Brown } 466c306944fSJed Brown 467d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetDim(DMPolytopeType ct) 468d71ae5a4SJacob Faibussowitsch { 469412e9a14SMatthew G. Knepley switch (ct) { 470d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT: 471d71ae5a4SJacob Faibussowitsch return 0; 472412e9a14SMatthew G. Knepley case DM_POLYTOPE_SEGMENT: 473d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT_PRISM_TENSOR: 474d71ae5a4SJacob Faibussowitsch return 1; 475412e9a14SMatthew G. Knepley case DM_POLYTOPE_TRIANGLE: 476412e9a14SMatthew G. Knepley case DM_POLYTOPE_QUADRILATERAL: 477d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_SEG_PRISM_TENSOR: 478476787b7SMatthew G. Knepley case DM_POLYTOPE_UNKNOWN_FACE: 479d71ae5a4SJacob Faibussowitsch return 2; 480412e9a14SMatthew G. Knepley case DM_POLYTOPE_TETRAHEDRON: 481412e9a14SMatthew G. Knepley case DM_POLYTOPE_HEXAHEDRON: 482412e9a14SMatthew G. Knepley case DM_POLYTOPE_TRI_PRISM: 483412e9a14SMatthew G. Knepley case DM_POLYTOPE_TRI_PRISM_TENSOR: 484412e9a14SMatthew G. Knepley case DM_POLYTOPE_QUAD_PRISM_TENSOR: 485d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_PYRAMID: 486476787b7SMatthew G. Knepley case DM_POLYTOPE_UNKNOWN_CELL: 487d71ae5a4SJacob Faibussowitsch return 3; 488d71ae5a4SJacob Faibussowitsch default: 489d71ae5a4SJacob Faibussowitsch return -1; 490412e9a14SMatthew G. Knepley } 491412e9a14SMatthew G. Knepley } 492412e9a14SMatthew G. Knepley 493d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetConeSize(DMPolytopeType ct) 494d71ae5a4SJacob Faibussowitsch { 495412e9a14SMatthew G. Knepley switch (ct) { 496d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT: 497d71ae5a4SJacob Faibussowitsch return 0; 498d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_SEGMENT: 499d71ae5a4SJacob Faibussowitsch return 2; 500d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT_PRISM_TENSOR: 501d71ae5a4SJacob Faibussowitsch return 2; 502d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRIANGLE: 503d71ae5a4SJacob Faibussowitsch return 3; 504d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUADRILATERAL: 505d71ae5a4SJacob Faibussowitsch return 4; 506d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_SEG_PRISM_TENSOR: 507d71ae5a4SJacob Faibussowitsch return 4; 508d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TETRAHEDRON: 509d71ae5a4SJacob Faibussowitsch return 4; 510d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_HEXAHEDRON: 511d71ae5a4SJacob Faibussowitsch return 6; 512d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRI_PRISM: 513d71ae5a4SJacob Faibussowitsch return 5; 514d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRI_PRISM_TENSOR: 515d71ae5a4SJacob Faibussowitsch return 5; 516d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUAD_PRISM_TENSOR: 517d71ae5a4SJacob Faibussowitsch return 6; 518d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_PYRAMID: 519d71ae5a4SJacob Faibussowitsch return 5; 520d71ae5a4SJacob Faibussowitsch default: 521d71ae5a4SJacob Faibussowitsch return -1; 522412e9a14SMatthew G. Knepley } 523412e9a14SMatthew G. Knepley } 524412e9a14SMatthew G. Knepley 525d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeGetNumVertices(DMPolytopeType ct) 526d71ae5a4SJacob Faibussowitsch { 527412e9a14SMatthew G. Knepley switch (ct) { 528d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT: 529d71ae5a4SJacob Faibussowitsch return 1; 530d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_SEGMENT: 531d71ae5a4SJacob Faibussowitsch return 2; 532d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT_PRISM_TENSOR: 533d71ae5a4SJacob Faibussowitsch return 2; 534d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRIANGLE: 535d71ae5a4SJacob Faibussowitsch return 3; 536d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUADRILATERAL: 537d71ae5a4SJacob Faibussowitsch return 4; 538d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_SEG_PRISM_TENSOR: 539d71ae5a4SJacob Faibussowitsch return 4; 540d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TETRAHEDRON: 541d71ae5a4SJacob Faibussowitsch return 4; 542d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_HEXAHEDRON: 543d71ae5a4SJacob Faibussowitsch return 8; 544d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRI_PRISM: 545d71ae5a4SJacob Faibussowitsch return 6; 546d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRI_PRISM_TENSOR: 547d71ae5a4SJacob Faibussowitsch return 6; 548d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUAD_PRISM_TENSOR: 549d71ae5a4SJacob Faibussowitsch return 8; 550d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_PYRAMID: 551d71ae5a4SJacob Faibussowitsch return 5; 552d71ae5a4SJacob Faibussowitsch default: 553d71ae5a4SJacob Faibussowitsch return -1; 554412e9a14SMatthew G. Knepley } 555412e9a14SMatthew G. Knepley } 556412e9a14SMatthew G. Knepley 557d71ae5a4SJacob Faibussowitsch static inline DMPolytopeType DMPolytopeTypeSimpleShape(PetscInt dim, PetscBool simplex) 558d71ae5a4SJacob Faibussowitsch { 5599371c9d4SSatish 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))); 5609318fe57SMatthew G. Knepley } 5619318fe57SMatthew G. Knepley 56285036b15SMatthew G. Knepley static inline PetscInt DMPolytopeTypeGetNumArrangements(DMPolytopeType ct) 563d71ae5a4SJacob Faibussowitsch { 564b5a892a1SMatthew G. Knepley switch (ct) { 565d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT: 566d71ae5a4SJacob Faibussowitsch return 1; 567d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_SEGMENT: 568d71ae5a4SJacob Faibussowitsch return 2; 569d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT_PRISM_TENSOR: 570d71ae5a4SJacob Faibussowitsch return 2; 571d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRIANGLE: 572d71ae5a4SJacob Faibussowitsch return 6; 573d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUADRILATERAL: 574d71ae5a4SJacob Faibussowitsch return 8; 575d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_SEG_PRISM_TENSOR: 576d71ae5a4SJacob Faibussowitsch return 4; 577d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TETRAHEDRON: 578d71ae5a4SJacob Faibussowitsch return 24; 579d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_HEXAHEDRON: 580d71ae5a4SJacob Faibussowitsch return 48; 581d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRI_PRISM: 582d71ae5a4SJacob Faibussowitsch return 12; 583d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRI_PRISM_TENSOR: 584d71ae5a4SJacob Faibussowitsch return 12; 585d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUAD_PRISM_TENSOR: 586d71ae5a4SJacob Faibussowitsch return 16; 587d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_PYRAMID: 588d71ae5a4SJacob Faibussowitsch return 8; 589d71ae5a4SJacob Faibussowitsch default: 590d71ae5a4SJacob Faibussowitsch return -1; 591b5a892a1SMatthew G. Knepley } 592b5a892a1SMatthew G. Knepley } 593b5a892a1SMatthew G. Knepley 594b5a892a1SMatthew G. Knepley /* An arrangement is a face order combined with an orientation for each face */ 59585036b15SMatthew G. Knepley static inline const PetscInt *DMPolytopeTypeGetArrangement(DMPolytopeType ct, PetscInt o) 596d71ae5a4SJacob Faibussowitsch { 597ef8b56bfSJed Brown static const PetscInt pntArr[1 * 2] = {0, 0}; 598b5a892a1SMatthew G. Knepley /* a: swap */ 5999371c9d4SSatish Balay static const PetscInt segArr[2 * 2 * 2] = {1, 0, 0, 0, /* -1: a */ 6009371c9d4SSatish Balay 0, 0, 1, 0, 6019371c9d4SSatish Balay /* 0: e */}; 602b5a892a1SMatthew G. Knepley /* a: swap first two 603b5a892a1SMatthew G. Knepley b: swap last two */ 6049371c9d4SSatish Balay static const PetscInt triArr[6 * 3 * 2] = {0, -1, 2, -1, 1, -1, /* -3: b */ 605b5a892a1SMatthew G. Knepley 2, -1, 1, -1, 0, -1, /* -2: aba */ 606b5a892a1SMatthew G. Knepley 1, -1, 0, -1, 2, -1, /* -1: a */ 607b5a892a1SMatthew G. Knepley 0, 0, 1, 0, 2, 0, /* 0: identity */ 608b5a892a1SMatthew G. Knepley 1, 0, 2, 0, 0, 0, /* 1: ba */ 6099371c9d4SSatish Balay 2, 0, 0, 0, 1, 0, 6109371c9d4SSatish Balay /* 2: ab */}; 611b5a892a1SMatthew G. Knepley /* a: forward cyclic permutation 612b5a892a1SMatthew G. Knepley b: swap first and last pairs */ 6139371c9d4SSatish Balay static const PetscInt quadArr[8 * 4 * 2] = {1, -1, 0, -1, 3, -1, 2, -1, /* -4: b */ 614b5a892a1SMatthew G. Knepley 0, -1, 3, -1, 2, -1, 1, -1, /* -3: b a^3 = a b */ 615b5a892a1SMatthew G. Knepley 3, -1, 2, -1, 1, -1, 0, -1, /* -2: b a^2 = a^2 b */ 616b5a892a1SMatthew G. Knepley 2, -1, 1, -1, 0, -1, 3, -1, /* -1: b a = a^3 b */ 617b5a892a1SMatthew G. Knepley 0, 0, 1, 0, 2, 0, 3, 0, /* 0: identity */ 618b5a892a1SMatthew G. Knepley 1, 0, 2, 0, 3, 0, 0, 0, /* 1: a */ 619b5a892a1SMatthew G. Knepley 2, 0, 3, 0, 0, 0, 1, 0, /* 2: a^2 */ 6209371c9d4SSatish Balay 3, 0, 0, 0, 1, 0, 2, 0, 6219371c9d4SSatish Balay /* 3: a^3 */}; 622b5a892a1SMatthew G. Knepley /* r: rotate 180 623b5a892a1SMatthew G. Knepley b: swap top and bottom segments */ 6249371c9d4SSatish Balay static const PetscInt tsegArr[4 * 4 * 2] = {1, -1, 0, -1, 3, -1, 2, -1, /* -2: r b */ 625b5a892a1SMatthew G. Knepley 0, -1, 1, -1, 3, 0, 2, 0, /* -1: r */ 626b5a892a1SMatthew G. Knepley 0, 0, 1, 0, 2, 0, 3, 0, /* 0: identity */ 6279371c9d4SSatish Balay 1, 0, 0, 0, 2, -1, 3, -1, 6289371c9d4SSatish Balay /* 1: b */}; 629b5a892a1SMatthew G. Knepley /* https://en.wikiversity.org/wiki/Symmetric_group_S4 */ 6309371c9d4SSatish Balay static const PetscInt tetArr[24 * 4 * 2] = {3, -2, 2, -3, 0, -1, 1, -1, /* -12: (1324) p22 */ 631b5a892a1SMatthew G. Knepley 3, -1, 1, -3, 2, -1, 0, -1, /* -11: (14) p21 */ 632b5a892a1SMatthew G. Knepley 3, -3, 0, -3, 1, -1, 2, -1, /* -10: (1234) p18 */ 633b5a892a1SMatthew G. Knepley 2, -1, 3, -1, 1, -3, 0, -2, /* -9: (1423) p17 */ 634b5a892a1SMatthew G. Knepley 2, -3, 0, -1, 3, -2, 1, -3, /* -8: (1342) p13 */ 635b5a892a1SMatthew G. Knepley 2, -2, 1, -2, 0, -2, 3, -2, /* -7: (24) p14 */ 636b5a892a1SMatthew G. Knepley 1, -2, 0, -2, 2, -2, 3, -1, /* -6: (34) p6 */ 637b5a892a1SMatthew G. Knepley 1, -1, 3, -3, 0, -3, 2, -2, /* -5: (1243) p10 */ 638b5a892a1SMatthew G. Knepley 1, -3, 2, -1, 3, -1, 0, -3, /* -4: (1432) p9 */ 639b5a892a1SMatthew G. Knepley 0, -3, 1, -1, 3, -3, 2, -3, /* -3: (12) p1 */ 640b5a892a1SMatthew G. Knepley 0, -2, 2, -2, 1, -2, 3, -3, /* -2: (23) p2 */ 641b5a892a1SMatthew G. Knepley 0, -1, 3, -2, 2, -3, 1, -2, /* -1: (13) p5 */ 642b5a892a1SMatthew G. Knepley 0, 0, 1, 0, 2, 0, 3, 0, /* 0: () p0 */ 643b5a892a1SMatthew G. Knepley 0, 1, 3, 1, 1, 2, 2, 0, /* 1: (123) p4 */ 644b5a892a1SMatthew G. Knepley 0, 2, 2, 1, 3, 0, 1, 2, /* 2: (132) p3 */ 645b5a892a1SMatthew G. Knepley 1, 2, 0, 1, 3, 1, 2, 2, /* 3: (12)(34) p7 */ 646b5a892a1SMatthew G. Knepley 1, 0, 2, 0, 0, 0, 3, 1, /* 4: (243) p8 */ 647b5a892a1SMatthew G. Knepley 1, 1, 3, 2, 2, 2, 0, 0, /* 5: (143) p11 */ 648b5a892a1SMatthew G. Knepley 2, 1, 3, 0, 0, 2, 1, 0, /* 6: (13)(24) p16 */ 649b5a892a1SMatthew G. Knepley 2, 2, 1, 1, 3, 2, 0, 2, /* 7: (142) p15 */ 650b5a892a1SMatthew G. Knepley 2, 0, 0, 0, 1, 0, 3, 2, /* 8: (234) p12 */ 651b5a892a1SMatthew G. Knepley 3, 2, 2, 2, 1, 1, 0, 1, /* 9: (14)(23) p23 */ 652b5a892a1SMatthew G. Knepley 3, 0, 0, 2, 2, 1, 1, 1, /* 10: (134) p19 */ 653b5a892a1SMatthew G. Knepley 3, 1, 1, 2, 0, 1, 2, 1 /* 11: (124) p20 */}; 654b5a892a1SMatthew G. Knepley /* Each rotation determines a permutation of the four diagonals, and this defines the isomorphism with S_4 */ 655ef8b56bfSJed Brown static const PetscInt hexArr[48 * 6 * 2] = { 656b5a892a1SMatthew G. Knepley 2, -3, 3, -2, 4, -2, 5, -3, 1, -3, 0, -1, /* -24: reflect bottom and use -3 on top */ 657b5a892a1SMatthew G. Knepley 4, -2, 5, -2, 0, -1, 1, -4, 3, -2, 2, -3, /* -23: reflect bottom and use -3 on top */ 658b5a892a1SMatthew G. Knepley 5, -3, 4, -1, 1, -2, 0, -3, 3, -4, 2, -1, /* -22: reflect bottom and use -3 on top */ 659b5a892a1SMatthew G. Knepley 3, -1, 2, -4, 4, -4, 5, -1, 0, -4, 1, -4, /* -21: reflect bottom and use -3 on top */ 660b5a892a1SMatthew G. Knepley 3, -3, 2, -2, 5, -1, 4, -4, 1, -1, 0, -3, /* -20: reflect bottom and use -3 on top */ 661b5a892a1SMatthew G. Knepley 4, -4, 5, -4, 1, -4, 0, -1, 2, -4, 3, -1, /* -19: reflect bottom and use -3 on top */ 662b5a892a1SMatthew G. Knepley 2, -1, 3, -4, 5, -3, 4, -2, 0, -2, 1, -2, /* -18: reflect bottom and use -3 on top */ 663b5a892a1SMatthew G. Knepley 5, -1, 4, -3, 0, -3, 1, -2, 2, -2, 3, -3, /* -17: reflect bottom and use -3 on top */ 664b5a892a1SMatthew G. Knepley 4, -3, 5, -1, 3, -2, 2, -4, 1, -4, 0, -4, /* -16: reflect bottom and use -3 on top */ 665b5a892a1SMatthew G. Knepley 5, -4, 4, -4, 3, -4, 2, -2, 0, -3, 1, -1, /* -15: reflect bottom and use -3 on top */ 666b5a892a1SMatthew G. Knepley 3, -4, 2, -1, 1, -1, 0, -4, 4, -4, 5, -4, /* -14: reflect bottom and use -3 on top */ 667b5a892a1SMatthew G. Knepley 2, -2, 3, -3, 0, -2, 1, -3, 4, -2, 5, -2, /* -13: reflect bottom and use -3 on top */ 668b5a892a1SMatthew G. Knepley 1, -3, 0, -1, 4, -1, 5, -4, 3, -1, 2, -4, /* -12: reflect bottom and use -3 on top */ 669b5a892a1SMatthew G. Knepley 1, -1, 0, -3, 5, -4, 4, -1, 2, -1, 3, -4, /* -11: reflect bottom and use -3 on top */ 670b5a892a1SMatthew G. Knepley 5, -2, 4, -2, 2, -2, 3, -4, 1, -2, 0, -2, /* -10: reflect bottom and use -3 on top */ 671b5a892a1SMatthew G. Knepley 1, -2, 0, -2, 2, -1, 3, -1, 4, -1, 5, -3, /* -9: reflect bottom and use -3 on top */ 672b5a892a1SMatthew G. Knepley 4, -1, 5, -3, 2, -4, 3, -2, 0, -1, 1, -3, /* -8: reflect bottom and use -3 on top */ 673b5a892a1SMatthew G. Knepley 3, -2, 2, -3, 0, -4, 1, -1, 5, -1, 4, -3, /* -7: reflect bottom and use -3 on top */ 674b5a892a1SMatthew G. Knepley 1, -4, 0, -4, 3, -1, 2, -1, 5, -4, 4, -4, /* -6: reflect bottom and use -3 on top */ 675b5a892a1SMatthew G. Knepley 2, -4, 3, -1, 1, -3, 0, -2, 5, -3, 4, -1, /* -5: reflect bottom and use -3 on top */ 676b5a892a1SMatthew G. Knepley 0, -4, 1, -4, 4, -3, 5, -2, 2, -3, 3, -2, /* -4: reflect bottom and use -3 on top */ 677b5a892a1SMatthew G. Knepley 0, -3, 1, -1, 3, -3, 2, -3, 4, -3, 5, -1, /* -3: reflect bottom and use -3 on top */ 678b5a892a1SMatthew G. Knepley 0, -2, 1, -2, 5, -2, 4, -3, 3, -3, 2, -2, /* -2: reflect bottom and use -3 on top */ 679b5a892a1SMatthew G. Knepley 0, -1, 1, -3, 2, -3, 3, -3, 5, -2, 4, -2, /* -1: reflect bottom and use -3 on top */ 680b5a892a1SMatthew G. Knepley 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, 5, 0, /* 0: identity */ 681b5a892a1SMatthew G. Knepley 0, 1, 1, 3, 5, 3, 4, 0, 2, 0, 3, 1, /* 1: 90 rotation about z */ 682b5a892a1SMatthew G. Knepley 0, 2, 1, 2, 3, 0, 2, 0, 5, 3, 4, 1, /* 2: 180 rotation about z */ 683b5a892a1SMatthew G. Knepley 0, 3, 1, 1, 4, 0, 5, 3, 3, 0, 2, 1, /* 3: 270 rotation about z */ 684b5a892a1SMatthew G. Knepley 2, 3, 3, 2, 1, 0, 0, 3, 4, 3, 5, 1, /* 4: 90 rotation about x */ 685b5a892a1SMatthew G. Knepley 1, 3, 0, 1, 3, 2, 2, 2, 4, 2, 5, 2, /* 5: 180 rotation about x */ 686b5a892a1SMatthew G. Knepley 3, 1, 2, 0, 0, 1, 1, 2, 4, 1, 5, 3, /* 6: 270 rotation about x */ 687b5a892a1SMatthew G. Knepley 4, 0, 5, 0, 2, 1, 3, 3, 1, 1, 0, 3, /* 7: 90 rotation about y */ 688b5a892a1SMatthew G. Knepley 1, 1, 0, 3, 2, 2, 3, 2, 5, 1, 4, 3, /* 8: 180 rotation about y */ 689b5a892a1SMatthew G. Knepley 5, 1, 4, 3, 2, 3, 3, 1, 0, 0, 1, 0, /* 9: 270 rotation about y */ 690b5a892a1SMatthew G. Knepley 1, 0, 0, 0, 5, 1, 4, 2, 3, 2, 2, 3, /* 10: 180 rotation about x+y */ 691b5a892a1SMatthew G. Knepley 1, 2, 0, 2, 4, 2, 5, 1, 2, 2, 3, 3, /* 11: 180 rotation about x-y */ 692b5a892a1SMatthew G. Knepley 2, 1, 3, 0, 0, 3, 1, 0, 5, 0, 4, 0, /* 12: 180 rotation about y+z */ 693b5a892a1SMatthew G. Knepley 3, 3, 2, 2, 1, 2, 0, 1, 5, 2, 4, 2, /* 13: 180 rotation about y-z */ 694b5a892a1SMatthew G. Knepley 5, 3, 4, 1, 3, 1, 2, 3, 1, 3, 0, 1, /* 14: 180 rotation about z+x */ 695b5a892a1SMatthew G. Knepley 4, 2, 5, 2, 3, 3, 2, 1, 0, 2, 1, 2, /* 15: 180 rotation about z-x */ 696b5a892a1SMatthew G. Knepley 5, 0, 4, 0, 0, 0, 1, 3, 3, 1, 2, 0, /* 16: 120 rotation about x+y+z (v0v6) */ 697b5a892a1SMatthew G. Knepley 2, 0, 3, 1, 5, 0, 4, 3, 1, 0, 0, 0, /* 17: 240 rotation about x+y+z (v0v6) */ 698b5a892a1SMatthew G. Knepley 4, 3, 5, 1, 1, 1, 0, 2, 3, 3, 2, 2, /* 18: 120 rotation about x+y-z (v4v2) */ 699b5a892a1SMatthew G. Knepley 3, 2, 2, 3, 5, 2, 4, 1, 0, 1, 1, 3, /* 19: 240 rotation about x+y-z (v4v2) */ 700b5a892a1SMatthew G. Knepley 3, 0, 2, 1, 4, 1, 5, 2, 1, 2, 0, 2, /* 20: 120 rotation about x-y+z (v1v5) */ 701b5a892a1SMatthew G. Knepley 5, 2, 4, 2, 1, 3, 0, 0, 2, 3, 3, 2, /* 21: 240 rotation about x-y+z (v1v5) */ 702b5a892a1SMatthew G. Knepley 4, 1, 5, 3, 0, 2, 1, 1, 2, 1, 3, 0, /* 22: 120 rotation about x-y-z (v7v3) */ 703b5a892a1SMatthew G. Knepley 2, 2, 3, 3, 4, 3, 5, 0, 0, 3, 1, 1, /* 23: 240 rotation about x-y-z (v7v3) */ 704b5a892a1SMatthew G. Knepley }; 705ef8b56bfSJed Brown static const PetscInt tripArr[12 * 5 * 2] = { 706b5a892a1SMatthew G. Knepley 1, -3, 0, -1, 3, -1, 4, -1, 2, -1, /* -6: reflect bottom and top */ 707b5a892a1SMatthew G. Knepley 1, -1, 0, -3, 4, -1, 2, -1, 3, -1, /* -5: reflect bottom and top */ 708b5a892a1SMatthew G. Knepley 1, -2, 0, -2, 2, -1, 3, -1, 4, -1, /* -4: reflect bottom and top */ 709b5a892a1SMatthew G. Knepley 0, -3, 1, -1, 3, -3, 2, -3, 4, -3, /* -3: reflect bottom and top */ 710b5a892a1SMatthew G. Knepley 0, -2, 1, -2, 4, -3, 3, -3, 2, -3, /* -2: reflect bottom and top */ 711b5a892a1SMatthew G. Knepley 0, -1, 1, -3, 2, -3, 4, -3, 3, -3, /* -1: reflect bottom and top */ 712b5a892a1SMatthew G. Knepley 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, /* 0: identity */ 713b5a892a1SMatthew G. Knepley 0, 1, 1, 2, 4, 0, 2, 0, 3, 0, /* 1: 120 rotation about z */ 714b5a892a1SMatthew G. Knepley 0, 2, 1, 1, 3, 0, 4, 0, 2, 0, /* 2: 240 rotation about z */ 715b5a892a1SMatthew G. Knepley 1, 1, 0, 2, 2, 2, 4, 2, 3, 2, /* 3: 180 rotation about y of 0 */ 716b5a892a1SMatthew G. Knepley 1, 0, 0, 0, 4, 2, 3, 2, 2, 2, /* 4: 180 rotation about y of 1 */ 717b5a892a1SMatthew G. Knepley 1, 2, 0, 1, 3, 2, 2, 2, 4, 2, /* 5: 180 rotation about y of 2 */ 718b5a892a1SMatthew G. Knepley }; 719b5a892a1SMatthew G. Knepley /* a: rotate 120 about z 720b5a892a1SMatthew G. Knepley b: swap top and bottom segments 721b5a892a1SMatthew G. Knepley r: reflect */ 722ef8b56bfSJed Brown static const PetscInt ttriArr[12 * 5 * 2] = { 723b5a892a1SMatthew G. Knepley 1, -3, 0, -3, 2, -2, 4, -2, 3, -2, /* -6: r b a^2 */ 724b5a892a1SMatthew G. Knepley 1, -2, 0, -2, 4, -2, 3, -2, 2, -2, /* -5: r b a */ 725b5a892a1SMatthew G. Knepley 1, -1, 0, -1, 3, -2, 2, -2, 4, -2, /* -4: r b */ 726b5a892a1SMatthew G. Knepley 0, -3, 1, -3, 2, -1, 4, -1, 3, -1, /* -3: r a^2 */ 727b5a892a1SMatthew G. Knepley 0, -2, 1, -2, 4, -1, 3, -1, 2, -1, /* -2: r a */ 728b5a892a1SMatthew G. Knepley 0, -1, 1, -1, 3, -1, 2, -1, 4, -1, /* -1: r */ 729b5a892a1SMatthew G. Knepley 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, /* 0: identity */ 730b5a892a1SMatthew G. Knepley 0, 1, 1, 1, 3, 0, 4, 0, 2, 0, /* 1: a */ 731b5a892a1SMatthew G. Knepley 0, 2, 1, 2, 4, 0, 2, 0, 3, 0, /* 2: a^2 */ 732b5a892a1SMatthew G. Knepley 1, 0, 0, 0, 2, 1, 3, 1, 4, 1, /* 3: b */ 733b5a892a1SMatthew G. Knepley 1, 1, 0, 1, 3, 1, 4, 1, 2, 1, /* 4: b a */ 734b5a892a1SMatthew G. Knepley 1, 2, 0, 2, 4, 1, 2, 1, 3, 1, /* 5: b a^2 */ 735b5a892a1SMatthew G. Knepley }; 736b5a892a1SMatthew G. Knepley /* a: rotate 90 about z 737b5a892a1SMatthew G. Knepley b: swap top and bottom segments 738b5a892a1SMatthew G. Knepley r: reflect */ 739ef8b56bfSJed Brown static const PetscInt tquadArr[16 * 6 * 2] = { 740b5a892a1SMatthew G. Knepley 1, -4, 0, -4, 3, -2, 2, -2, 5, -2, 4, -2, /* -8: r b a^3 */ 741b5a892a1SMatthew G. Knepley 1, -3, 0, -3, 2, -2, 5, -2, 4, -2, 3, -2, /* -7: r b a^2 */ 742b5a892a1SMatthew G. Knepley 1, -2, 0, -2, 5, -2, 4, -2, 3, -2, 2, -2, /* -6: r b a */ 743b5a892a1SMatthew G. Knepley 1, -1, 0, -1, 4, -2, 3, -2, 2, -2, 5, -2, /* -5: r b */ 744b5a892a1SMatthew G. Knepley 0, -4, 1, -4, 3, -1, 2, -1, 5, -1, 4, -1, /* -4: r a^3 */ 745b5a892a1SMatthew G. Knepley 0, -3, 1, -3, 2, -1, 5, -1, 4, -1, 3, -1, /* -3: r a^2 */ 746b5a892a1SMatthew G. Knepley 0, -2, 1, -2, 5, -1, 4, -1, 3, -1, 2, -1, /* -2: r a */ 747b5a892a1SMatthew G. Knepley 0, -1, 1, -1, 4, -1, 3, -1, 2, -1, 5, -1, /* -1: r */ 748b5a892a1SMatthew G. Knepley 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, 5, 0, /* 0: identity */ 749b5a892a1SMatthew G. Knepley 0, 1, 1, 1, 3, 0, 4, 0, 5, 0, 2, 0, /* 1: a */ 750b5a892a1SMatthew G. Knepley 0, 2, 1, 2, 4, 0, 5, 0, 2, 0, 3, 0, /* 2: a^2 */ 751b5a892a1SMatthew G. Knepley 0, 3, 1, 3, 5, 0, 2, 0, 3, 0, 4, 0, /* 3: a^3 */ 752b5a892a1SMatthew G. Knepley 1, 0, 0, 0, 2, 1, 3, 1, 4, 1, 5, 1, /* 4: b */ 753b5a892a1SMatthew G. Knepley 1, 1, 0, 1, 3, 1, 4, 1, 5, 1, 2, 1, /* 5: b a */ 754b5a892a1SMatthew G. Knepley 1, 2, 0, 2, 4, 1, 5, 1, 2, 1, 3, 1, /* 6: b a^2 */ 755b5a892a1SMatthew G. Knepley 1, 3, 0, 3, 5, 1, 2, 1, 3, 1, 4, 1, /* 7: b a^3 */ 756b5a892a1SMatthew G. Knepley }; 757ef8b56bfSJed Brown static const PetscInt pyrArr[8 * 5 * 2] = { 758b5a892a1SMatthew G. Knepley 0, -4, 2, -3, 1, -3, 4, -3, 3, -3, /* -4: Reflect bottom face */ 759b5a892a1SMatthew G. Knepley 0, -3, 3, -3, 2, -3, 1, -3, 4, -3, /* -3: Reflect bottom face */ 760b5a892a1SMatthew G. Knepley 0, -2, 4, -3, 3, -3, 2, -3, 1, -3, /* -2: Reflect bottom face */ 761b5a892a1SMatthew G. Knepley 0, -1, 1, -3, 4, -3, 3, -3, 2, -3, /* -1: Reflect bottom face */ 762b5a892a1SMatthew G. Knepley 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, /* 0: identity */ 763b5a892a1SMatthew G. Knepley 0, 1, 4, 0, 1, 0, 2, 0, 3, 0, /* 1: 90 rotation about z */ 764b5a892a1SMatthew G. Knepley 0, 2, 3, 0, 4, 0, 1, 0, 2, 0, /* 2: 180 rotation about z */ 765b5a892a1SMatthew G. Knepley 0, 3, 2, 0, 3, 0, 4, 0, 1, 0, /* 3: 270 rotation about z */ 766b5a892a1SMatthew G. Knepley }; 767b5a892a1SMatthew G. Knepley switch (ct) { 768d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT: 769d71ae5a4SJacob Faibussowitsch return pntArr; 770d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_SEGMENT: 771d71ae5a4SJacob Faibussowitsch return &segArr[(o + 1) * 2 * 2]; 772d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT_PRISM_TENSOR: 773d71ae5a4SJacob Faibussowitsch return &segArr[(o + 1) * 2 * 2]; 774d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRIANGLE: 775d71ae5a4SJacob Faibussowitsch return &triArr[(o + 3) * 3 * 2]; 776d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUADRILATERAL: 777d71ae5a4SJacob Faibussowitsch return &quadArr[(o + 4) * 4 * 2]; 778d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_SEG_PRISM_TENSOR: 779d71ae5a4SJacob Faibussowitsch return &tsegArr[(o + 2) * 4 * 2]; 780d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TETRAHEDRON: 781d71ae5a4SJacob Faibussowitsch return &tetArr[(o + 12) * 4 * 2]; 782d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_HEXAHEDRON: 783d71ae5a4SJacob Faibussowitsch return &hexArr[(o + 24) * 6 * 2]; 784d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRI_PRISM: 785d71ae5a4SJacob Faibussowitsch return &tripArr[(o + 6) * 5 * 2]; 786d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRI_PRISM_TENSOR: 787d71ae5a4SJacob Faibussowitsch return &ttriArr[(o + 6) * 5 * 2]; 788d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUAD_PRISM_TENSOR: 789d71ae5a4SJacob Faibussowitsch return &tquadArr[(o + 8) * 6 * 2]; 790d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_PYRAMID: 791d71ae5a4SJacob Faibussowitsch return &pyrArr[(o + 4) * 5 * 2]; 792d71ae5a4SJacob Faibussowitsch default: 793f22e26b7SPierre Jolivet return PETSC_NULLPTR; 794b5a892a1SMatthew G. Knepley } 795b5a892a1SMatthew G. Knepley } 796b5a892a1SMatthew G. Knepley 797d5b43468SJose E. Roman /* A vertex arrangement is a vertex order */ 79885036b15SMatthew G. Knepley static inline const PetscInt *DMPolytopeTypeGetVertexArrangement(DMPolytopeType ct, PetscInt o) 799d71ae5a4SJacob Faibussowitsch { 800ef8b56bfSJed Brown static const PetscInt pntVerts[1] = {0}; 8019371c9d4SSatish Balay static const PetscInt segVerts[2 * 2] = {1, 0, 0, 1}; 8029371c9d4SSatish 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}; 8039371c9d4SSatish 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}; 8049371c9d4SSatish Balay static const PetscInt tsegVerts[4 * 4] = {3, 2, 1, 0, 1, 0, 3, 2, 0, 1, 2, 3, 2, 3, 0, 1}; 8059371c9d4SSatish Balay static const PetscInt tetVerts[24 * 4] = {2, 3, 1, 0, /* -12: (1324) p22 */ 806b5a892a1SMatthew G. Knepley 3, 1, 2, 0, /* -11: (14) p21 */ 807b5a892a1SMatthew G. Knepley 1, 2, 3, 0, /* -10: (1234) p18 */ 808b5a892a1SMatthew G. Knepley 3, 2, 0, 1, /* -9: (1423) p17 */ 809b5a892a1SMatthew G. Knepley 2, 0, 3, 1, /* -8: (1342) p13 */ 810b5a892a1SMatthew G. Knepley 0, 3, 2, 1, /* -7: (24) p14 */ 811b5a892a1SMatthew G. Knepley 0, 1, 3, 2, /* -6: (34) p6 */ 812b5a892a1SMatthew G. Knepley 1, 3, 0, 2, /* -5: (1243) p10 */ 813b5a892a1SMatthew G. Knepley 3, 0, 1, 2, /* -4: (1432 p9 */ 814b5a892a1SMatthew G. Knepley 1, 0, 2, 3, /* -3: (12) p1 */ 815b5a892a1SMatthew G. Knepley 0, 2, 1, 3, /* -2: (23) p2 */ 816b5a892a1SMatthew G. Knepley 2, 1, 0, 3, /* -1: (13) p5 */ 817b5a892a1SMatthew G. Knepley 0, 1, 2, 3, /* 0: () p0 */ 818b5a892a1SMatthew G. Knepley 1, 2, 0, 3, /* 1: (123) p4 */ 819b5a892a1SMatthew G. Knepley 2, 0, 1, 3, /* 2: (132) p3 */ 820b5a892a1SMatthew G. Knepley 1, 0, 3, 2, /* 3: (12)(34) p7 */ 821b5a892a1SMatthew G. Knepley 0, 3, 1, 2, /* 4: (243) p8 */ 822b5a892a1SMatthew G. Knepley 3, 1, 0, 2, /* 5: (143) p11 */ 823b5a892a1SMatthew G. Knepley 2, 3, 0, 1, /* 6: (13)(24) p16 */ 824b5a892a1SMatthew G. Knepley 3, 0, 2, 1, /* 7: (142) p15 */ 825b5a892a1SMatthew G. Knepley 0, 2, 3, 1, /* 8: (234) p12 */ 826b5a892a1SMatthew G. Knepley 3, 2, 1, 0, /* 9: (14)(23) p23 */ 827b5a892a1SMatthew G. Knepley 2, 1, 3, 0, /* 10: (134) p19 */ 828b5a892a1SMatthew G. Knepley 1, 3, 2, 0 /* 11: (124) p20 */}; 829ef8b56bfSJed Brown static const PetscInt hexVerts[48 * 8] = { 830b5a892a1SMatthew G. Knepley 3, 0, 4, 5, 2, 6, 7, 1, /* -24: reflected 23 */ 831b5a892a1SMatthew G. Knepley 3, 5, 6, 2, 0, 1, 7, 4, /* -23: reflected 22 */ 832b5a892a1SMatthew G. Knepley 4, 0, 1, 7, 5, 6, 2, 3, /* -22: reflected 21 */ 833b5a892a1SMatthew G. Knepley 6, 7, 1, 2, 5, 3, 0, 4, /* -21: reflected 20 */ 834b5a892a1SMatthew G. Knepley 1, 2, 6, 7, 0, 4, 5, 3, /* -20: reflected 19 */ 835b5a892a1SMatthew G. Knepley 6, 2, 3, 5, 7, 4, 0, 1, /* -19: reflected 18 */ 836b5a892a1SMatthew G. Knepley 4, 5, 3, 0, 7, 1, 2, 6, /* -18: reflected 17 */ 837b5a892a1SMatthew G. Knepley 1, 7, 4, 0, 2, 3, 5, 6, /* -17: reflected 16 */ 838b5a892a1SMatthew G. Knepley 2, 3, 5, 6, 1, 7, 4, 0, /* -16: reflected 15 */ 839b5a892a1SMatthew G. Knepley 7, 4, 0, 1, 6, 2, 3, 5, /* -15: reflected 14 */ 840b5a892a1SMatthew G. Knepley 7, 1, 2, 6, 4, 5, 3, 0, /* -14: reflected 13 */ 841b5a892a1SMatthew G. Knepley 0, 4, 5, 3, 1, 2, 6, 7, /* -13: reflected 12 */ 842b5a892a1SMatthew G. Knepley 5, 4, 7, 6, 3, 2, 1, 0, /* -12: reflected 11 */ 843b5a892a1SMatthew G. Knepley 7, 6, 5, 4, 1, 0, 3, 2, /* -11: reflected 10 */ 844b5a892a1SMatthew G. Knepley 0, 1, 7, 4, 3, 5, 6, 2, /* -10: reflected 9 */ 845b5a892a1SMatthew G. Knepley 4, 7, 6, 5, 0, 3, 2, 1, /* -9: reflected 8 */ 846b5a892a1SMatthew G. Knepley 5, 6, 2, 3, 4, 0, 1, 7, /* -8: reflected 7 */ 847b5a892a1SMatthew G. Knepley 2, 6, 7, 1, 3, 0, 4, 5, /* -7: reflected 6 */ 848b5a892a1SMatthew G. Knepley 6, 5, 4, 7, 2, 1, 0, 3, /* -6: reflected 5 */ 849b5a892a1SMatthew G. Knepley 5, 3, 0, 4, 6, 7, 1, 2, /* -5: reflected 4 */ 850b5a892a1SMatthew G. Knepley 2, 1, 0, 3, 6, 5, 4, 7, /* -4: reflected 3 */ 851b5a892a1SMatthew G. Knepley 1, 0, 3, 2, 7, 6, 5, 4, /* -3: reflected 2 */ 852b5a892a1SMatthew G. Knepley 0, 3, 2, 1, 4, 7, 6, 5, /* -2: reflected 1 */ 853b5a892a1SMatthew G. Knepley 3, 2, 1, 0, 5, 4, 7, 6, /* -1: reflected 0 */ 854b5a892a1SMatthew G. Knepley 0, 1, 2, 3, 4, 5, 6, 7, /* 0: identity */ 855b5a892a1SMatthew G. Knepley 1, 2, 3, 0, 7, 4, 5, 6, /* 1: 90 rotation about z */ 856b5a892a1SMatthew G. Knepley 2, 3, 0, 1, 6, 7, 4, 5, /* 2: 180 rotation about z */ 857b5a892a1SMatthew G. Knepley 3, 0, 1, 2, 5, 6, 7, 4, /* 3: 270 rotation about z */ 858b5a892a1SMatthew G. Knepley 4, 0, 3, 5, 7, 6, 2, 1, /* 4: 90 rotation about x */ 859b5a892a1SMatthew G. Knepley 7, 4, 5, 6, 1, 2, 3, 0, /* 5: 180 rotation about x */ 860b5a892a1SMatthew G. Knepley 1, 7, 6, 2, 0, 3, 5, 4, /* 6: 270 rotation about x */ 861b5a892a1SMatthew G. Knepley 3, 2, 6, 5, 0, 4, 7, 1, /* 7: 90 rotation about y */ 862b5a892a1SMatthew G. Knepley 5, 6, 7, 4, 3, 0, 1, 2, /* 8: 180 rotation about y */ 863b5a892a1SMatthew G. Knepley 4, 7, 1, 0, 5, 3, 2, 6, /* 9: 270 rotation about y */ 864b5a892a1SMatthew G. Knepley 4, 5, 6, 7, 0, 1, 2, 3, /* 10: 180 rotation about x+y */ 865b5a892a1SMatthew G. Knepley 6, 7, 4, 5, 2, 3, 0, 1, /* 11: 180 rotation about x-y */ 866b5a892a1SMatthew G. Knepley 3, 5, 4, 0, 2, 1, 7, 6, /* 12: 180 rotation about y+z */ 867b5a892a1SMatthew G. Knepley 6, 2, 1, 7, 5, 4, 0, 3, /* 13: 180 rotation about y-z */ 868b5a892a1SMatthew G. Knepley 1, 0, 4, 7, 2, 6, 5, 3, /* 14: 180 rotation about z+x */ 869b5a892a1SMatthew G. Knepley 6, 5, 3, 2, 7, 1, 0, 4, /* 15: 180 rotation about z-x */ 870b5a892a1SMatthew G. Knepley 0, 4, 7, 1, 3, 2, 6, 5, /* 16: 120 rotation about x+y+z (v0v6) */ 871b5a892a1SMatthew G. Knepley 0, 3, 5, 4, 1, 7, 6, 2, /* 17: 240 rotation about x+y+z (v0v6) */ 872b5a892a1SMatthew G. Knepley 5, 3, 2, 6, 4, 7, 1, 0, /* 18: 120 rotation about x+y-z (v4v2) */ 873b5a892a1SMatthew G. Knepley 7, 6, 2, 1, 4, 0, 3, 5, /* 19: 240 rotation about x+y-z (v4v2) */ 874b5a892a1SMatthew G. Knepley 2, 1, 7, 6, 3, 5, 4, 0, /* 20: 120 rotation about x-y+z (v1v5) */ 875b5a892a1SMatthew G. Knepley 7, 1, 0, 4, 6, 5, 3, 2, /* 21: 240 rotation about x-y+z (v1v5) */ 876b5a892a1SMatthew G. Knepley 2, 6, 5, 3, 1, 0, 4, 7, /* 22: 120 rotation about x-y-z (v7v3) */ 877b5a892a1SMatthew G. Knepley 5, 4, 0, 3, 6, 2, 1, 7, /* 23: 240 rotation about x-y-z (v7v3) */ 878b5a892a1SMatthew G. Knepley }; 879ef8b56bfSJed Brown static const PetscInt tripVerts[12 * 6] = { 880b5a892a1SMatthew G. Knepley 4, 3, 5, 2, 1, 0, /* -6: reflect bottom and top */ 881b5a892a1SMatthew G. Knepley 5, 4, 3, 1, 0, 2, /* -5: reflect bottom and top */ 882b5a892a1SMatthew G. Knepley 3, 5, 4, 0, 2, 1, /* -4: reflect bottom and top */ 883b5a892a1SMatthew G. Knepley 1, 0, 2, 5, 4, 3, /* -3: reflect bottom and top */ 884b5a892a1SMatthew G. Knepley 0, 2, 1, 3, 5, 4, /* -2: reflect bottom and top */ 885b5a892a1SMatthew G. Knepley 2, 1, 0, 4, 3, 5, /* -1: reflect bottom and top */ 886b5a892a1SMatthew G. Knepley 0, 1, 2, 3, 4, 5, /* 0: identity */ 887b5a892a1SMatthew G. Knepley 1, 2, 0, 5, 3, 4, /* 1: 120 rotation about z */ 888b5a892a1SMatthew G. Knepley 2, 0, 1, 4, 5, 3, /* 2: 240 rotation about z */ 889b5a892a1SMatthew G. Knepley 4, 5, 3, 2, 0, 1, /* 3: 180 rotation about y of 0 */ 890b5a892a1SMatthew G. Knepley 3, 4, 5, 0, 1, 2, /* 4: 180 rotation about y of 1 */ 891b5a892a1SMatthew G. Knepley 5, 3, 4, 1, 2, 0, /* 5: 180 rotation about y of 2 */ 892b5a892a1SMatthew G. Knepley }; 893ef8b56bfSJed Brown static const PetscInt ttriVerts[12 * 6] = { 894b5a892a1SMatthew G. Knepley 4, 3, 5, 1, 0, 2, /* -6: r b a^2 */ 895b5a892a1SMatthew G. Knepley 3, 5, 4, 0, 2, 1, /* -5: r b a */ 896b5a892a1SMatthew G. Knepley 5, 4, 3, 2, 1, 0, /* -4: r b */ 897b5a892a1SMatthew G. Knepley 1, 0, 2, 4, 3, 5, /* -3: r a^2 */ 898b5a892a1SMatthew G. Knepley 0, 2, 1, 3, 5, 4, /* -2: r a */ 899b5a892a1SMatthew G. Knepley 2, 1, 0, 5, 4, 3, /* -1: r */ 900b5a892a1SMatthew G. Knepley 0, 1, 2, 3, 4, 5, /* 0: identity */ 901b5a892a1SMatthew G. Knepley 1, 2, 0, 4, 5, 3, /* 1: a */ 902b5a892a1SMatthew G. Knepley 2, 0, 1, 5, 3, 4, /* 2: a^2 */ 903b5a892a1SMatthew G. Knepley 3, 4, 5, 0, 1, 2, /* 3: b */ 904b5a892a1SMatthew G. Knepley 4, 5, 3, 1, 2, 0, /* 4: b a */ 905b5a892a1SMatthew G. Knepley 5, 3, 4, 2, 0, 1, /* 5: b a^2 */ 906b5a892a1SMatthew G. Knepley }; 907b5a892a1SMatthew G. Knepley /* a: rotate 90 about z 908b5a892a1SMatthew G. Knepley b: swap top and bottom segments 909b5a892a1SMatthew G. Knepley r: reflect */ 910ef8b56bfSJed Brown static const PetscInt tquadVerts[16 * 8] = { 911b5a892a1SMatthew G. Knepley 6, 5, 4, 7, 2, 1, 0, 3, /* -8: r b a^3 */ 912b5a892a1SMatthew G. Knepley 5, 4, 7, 6, 1, 0, 3, 2, /* -7: r b a^2 */ 913b5a892a1SMatthew G. Knepley 4, 7, 6, 5, 0, 3, 2, 1, /* -6: r b a */ 914b5a892a1SMatthew G. Knepley 7, 6, 5, 4, 3, 2, 1, 0, /* -5: r b */ 915b5a892a1SMatthew G. Knepley 2, 1, 0, 3, 6, 5, 4, 7, /* -4: r a^3 */ 916b5a892a1SMatthew G. Knepley 1, 0, 3, 2, 5, 4, 7, 6, /* -3: r a^2 */ 917b5a892a1SMatthew G. Knepley 0, 3, 2, 1, 4, 7, 6, 5, /* -2: r a */ 918b5a892a1SMatthew G. Knepley 3, 2, 1, 0, 7, 6, 5, 4, /* -1: r */ 919b5a892a1SMatthew G. Knepley 0, 1, 2, 3, 4, 5, 6, 7, /* 0: identity */ 920b5a892a1SMatthew G. Knepley 1, 2, 3, 0, 5, 6, 7, 4, /* 1: a */ 921b5a892a1SMatthew G. Knepley 2, 3, 0, 1, 6, 7, 4, 5, /* 2: a^2 */ 922b5a892a1SMatthew G. Knepley 3, 0, 1, 2, 7, 4, 5, 6, /* 3: a^3 */ 923b5a892a1SMatthew G. Knepley 4, 5, 6, 7, 0, 1, 2, 3, /* 4: b */ 924b5a892a1SMatthew G. Knepley 5, 6, 7, 4, 1, 2, 3, 0, /* 5: b a */ 925b5a892a1SMatthew G. Knepley 6, 7, 4, 5, 2, 3, 0, 1, /* 6: b a^2 */ 926b5a892a1SMatthew G. Knepley 7, 4, 5, 6, 3, 0, 1, 2, /* 7: b a^3 */ 927b5a892a1SMatthew G. Knepley }; 928ef8b56bfSJed Brown static const PetscInt pyrVerts[8 * 5] = { 929b5a892a1SMatthew G. Knepley 2, 1, 0, 3, 4, /* -4: Reflect bottom face */ 930b5a892a1SMatthew G. Knepley 1, 0, 3, 2, 4, /* -3: Reflect bottom face */ 931b5a892a1SMatthew G. Knepley 0, 3, 2, 1, 4, /* -2: Reflect bottom face */ 932b5a892a1SMatthew G. Knepley 3, 2, 1, 0, 4, /* -1: Reflect bottom face */ 933b5a892a1SMatthew G. Knepley 0, 1, 2, 3, 4, /* 0: identity */ 934b5a892a1SMatthew G. Knepley 1, 2, 3, 0, 4, /* 1: 90 rotation about z */ 935b5a892a1SMatthew G. Knepley 2, 3, 0, 1, 4, /* 2: 180 rotation about z */ 936b5a892a1SMatthew G. Knepley 3, 0, 1, 2, 4, /* 3: 270 rotation about z */ 937b5a892a1SMatthew G. Knepley }; 938b5a892a1SMatthew G. Knepley switch (ct) { 939d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT: 940d71ae5a4SJacob Faibussowitsch return pntVerts; 941d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_SEGMENT: 942d71ae5a4SJacob Faibussowitsch return &segVerts[(o + 1) * 2]; 943d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT_PRISM_TENSOR: 944d71ae5a4SJacob Faibussowitsch return &segVerts[(o + 1) * 2]; 945d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRIANGLE: 946d71ae5a4SJacob Faibussowitsch return &triVerts[(o + 3) * 3]; 947d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUADRILATERAL: 948d71ae5a4SJacob Faibussowitsch return &quadVerts[(o + 4) * 4]; 949d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_SEG_PRISM_TENSOR: 950d71ae5a4SJacob Faibussowitsch return &tsegVerts[(o + 2) * 4]; 951d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TETRAHEDRON: 952d71ae5a4SJacob Faibussowitsch return &tetVerts[(o + 12) * 4]; 953d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_HEXAHEDRON: 954d71ae5a4SJacob Faibussowitsch return &hexVerts[(o + 24) * 8]; 955d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRI_PRISM: 956d71ae5a4SJacob Faibussowitsch return &tripVerts[(o + 6) * 6]; 957d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRI_PRISM_TENSOR: 958d71ae5a4SJacob Faibussowitsch return &ttriVerts[(o + 6) * 6]; 959d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUAD_PRISM_TENSOR: 960d71ae5a4SJacob Faibussowitsch return &tquadVerts[(o + 8) * 8]; 961d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_PYRAMID: 962d71ae5a4SJacob Faibussowitsch return &pyrVerts[(o + 4) * 5]; 963d71ae5a4SJacob Faibussowitsch default: 964f22e26b7SPierre Jolivet return PETSC_NULLPTR; 965b5a892a1SMatthew G. Knepley } 966b5a892a1SMatthew G. Knepley } 967b5a892a1SMatthew G. Knepley 968b5a892a1SMatthew G. Knepley /* This is orientation o1 acting on orientation o2 */ 969d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeComposeOrientation(DMPolytopeType ct, PetscInt o1, PetscInt o2) 970d71ae5a4SJacob Faibussowitsch { 9719371c9d4SSatish Balay static const PetscInt segMult[2 * 2] = {0, -1, -1, 0}; 9729371c9d4SSatish 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}; 9739371c9d4SSatish 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, 9749371c9d4SSatish 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}; 9759371c9d4SSatish Balay static const PetscInt tsegMult[4 * 4] = {0, 1, -2, -1, 1, 0, -1, -2, -2, -1, 0, 1, -1, -2, 1, 0}; 976ef8b56bfSJed Brown static const PetscInt tetMult[24 * 24] = { 9779371c9d4SSatish 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, 9789371c9d4SSatish 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, 9799371c9d4SSatish 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, 9809371c9d4SSatish 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, 9819371c9d4SSatish 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, 9829371c9d4SSatish 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, 9839371c9d4SSatish 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, 9849371c9d4SSatish 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, 9859371c9d4SSatish 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, 9869371c9d4SSatish 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, 9879371c9d4SSatish 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, 9889371c9d4SSatish 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, 989b5a892a1SMatthew G. Knepley }; 990ef8b56bfSJed Brown static const PetscInt hexMult[48 * 48] = { 991b5a892a1SMatthew 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, 992b5a892a1SMatthew 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, 993b5a892a1SMatthew 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, 994b5a892a1SMatthew 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, 995b5a892a1SMatthew 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, 996b5a892a1SMatthew 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, 997b5a892a1SMatthew 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, 998b5a892a1SMatthew 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, 999b5a892a1SMatthew 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, 1000b5a892a1SMatthew 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, 1001b5a892a1SMatthew 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, 1002b5a892a1SMatthew 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, 1003b5a892a1SMatthew 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, 1004b5a892a1SMatthew 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, 1005b5a892a1SMatthew 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, 1006b5a892a1SMatthew 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, 1007b5a892a1SMatthew 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, 1008b5a892a1SMatthew 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, 1009b5a892a1SMatthew 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, 1010b5a892a1SMatthew 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, 1011b5a892a1SMatthew 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, 1012b5a892a1SMatthew 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, 1013b5a892a1SMatthew 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, 1014b5a892a1SMatthew 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, 1015b5a892a1SMatthew 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, 1016b5a892a1SMatthew 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, 1017b5a892a1SMatthew 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, 1018b5a892a1SMatthew 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, 1019b5a892a1SMatthew 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, 1020b5a892a1SMatthew 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, 1021b5a892a1SMatthew 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, 1022b5a892a1SMatthew 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, 1023b5a892a1SMatthew 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, 1024b5a892a1SMatthew 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, 1025b5a892a1SMatthew 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, 1026b5a892a1SMatthew 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, 1027b5a892a1SMatthew 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, 1028b5a892a1SMatthew 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, 1029b5a892a1SMatthew 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, 1030b5a892a1SMatthew 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, 1031b5a892a1SMatthew 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, 1032b5a892a1SMatthew 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, 1033b5a892a1SMatthew 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, 1034b5a892a1SMatthew 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, 1035b5a892a1SMatthew 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, 1036b5a892a1SMatthew 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, 1037b5a892a1SMatthew 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, 1038b5a892a1SMatthew 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, 1039b5a892a1SMatthew G. Knepley }; 1040ef8b56bfSJed Brown static const PetscInt tripMult[12 * 12] = { 10419371c9d4SSatish 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, 10429371c9d4SSatish 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, 10439371c9d4SSatish 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, 1044b5a892a1SMatthew G. Knepley }; 1045ef8b56bfSJed Brown static const PetscInt ttriMult[12 * 12] = { 10469371c9d4SSatish 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, 10479371c9d4SSatish 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, 10489371c9d4SSatish 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, 1049b5a892a1SMatthew G. Knepley }; 1050ef8b56bfSJed Brown static const PetscInt tquadMult[16 * 16] = { 10519371c9d4SSatish 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, 10529371c9d4SSatish 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, 10539371c9d4SSatish 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, 10549371c9d4SSatish 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, 10559371c9d4SSatish 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, 1056b5a892a1SMatthew G. Knepley }; 1057ef8b56bfSJed Brown static const PetscInt pyrMult[8 * 8] = { 10589371c9d4SSatish 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, 1059b5a892a1SMatthew G. Knepley }; 1060b5a892a1SMatthew G. Knepley switch (ct) { 1061d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT: 1062d71ae5a4SJacob Faibussowitsch return 0; 1063b5a892a1SMatthew G. Knepley case DM_POLYTOPE_SEGMENT: 1064d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT_PRISM_TENSOR: 1065d71ae5a4SJacob Faibussowitsch return segMult[(o1 + 1) * 2 + o2 + 1]; 1066d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRIANGLE: 1067d71ae5a4SJacob Faibussowitsch return triMult[(o1 + 3) * 6 + o2 + 3]; 1068d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUADRILATERAL: 1069d71ae5a4SJacob Faibussowitsch return quadMult[(o1 + 4) * 8 + o2 + 4]; 1070d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_SEG_PRISM_TENSOR: 1071d71ae5a4SJacob Faibussowitsch return tsegMult[(o1 + 2) * 4 + o2 + 2]; 1072d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TETRAHEDRON: 1073d71ae5a4SJacob Faibussowitsch return tetMult[(o1 + 12) * 24 + o2 + 12]; 1074d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_HEXAHEDRON: 1075d71ae5a4SJacob Faibussowitsch return hexMult[(o1 + 24) * 48 + o2 + 24]; 1076d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRI_PRISM: 1077d71ae5a4SJacob Faibussowitsch return tripMult[(o1 + 6) * 12 + o2 + 6]; 1078d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRI_PRISM_TENSOR: 1079d71ae5a4SJacob Faibussowitsch return ttriMult[(o1 + 6) * 12 + o2 + 6]; 1080d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUAD_PRISM_TENSOR: 1081d71ae5a4SJacob Faibussowitsch return tquadMult[(o1 + 8) * 16 + o2 + 8]; 1082d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_PYRAMID: 1083d71ae5a4SJacob Faibussowitsch return pyrMult[(o1 + 4) * 8 + o2 + 4]; 1084d71ae5a4SJacob Faibussowitsch default: 1085d71ae5a4SJacob Faibussowitsch return 0; 1086b5a892a1SMatthew G. Knepley } 1087b5a892a1SMatthew G. Knepley } 1088b5a892a1SMatthew G. Knepley 1089b5a892a1SMatthew G. Knepley /* This is orientation o1 acting on orientation o2^{-1} */ 1090d71ae5a4SJacob Faibussowitsch static inline PetscInt DMPolytopeTypeComposeOrientationInv(DMPolytopeType ct, PetscInt o1, PetscInt o2) 1091d71ae5a4SJacob Faibussowitsch { 1092ef8b56bfSJed Brown static const PetscInt triInv[6] = {-3, -2, -1, 0, 2, 1}; 1093ef8b56bfSJed Brown static const PetscInt quadInv[8] = {-4, -3, -2, -1, 0, 3, 2, 1}; 1094ef8b56bfSJed 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}; 10959371c9d4SSatish 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}; 1096ef8b56bfSJed Brown static const PetscInt tripInv[12] = {-5, -6, -4, -3, -2, -1, 0, 2, 1, 3, 4, 5}; 1097ef8b56bfSJed Brown static const PetscInt ttriInv[12] = {-6, -5, -4, -3, -2, -1, 0, 2, 1, 3, 5, 4}; 1098ef8b56bfSJed Brown static const PetscInt tquadInv[16] = {-8, -7, -6, -5, -4, -3, -2, -1, 0, 3, 2, 1, 4, 7, 6, 5}; 1099ef8b56bfSJed Brown static const PetscInt pyrInv[8] = {-4, -3, -2, -1, 0, 3, 2, 1}; 1100b5a892a1SMatthew G. Knepley switch (ct) { 1101d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT: 1102d71ae5a4SJacob Faibussowitsch return 0; 1103b5a892a1SMatthew G. Knepley case DM_POLYTOPE_SEGMENT: 1104d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_POINT_PRISM_TENSOR: 1105d71ae5a4SJacob Faibussowitsch return DMPolytopeTypeComposeOrientation(ct, o1, o2); 1106d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRIANGLE: 1107d71ae5a4SJacob Faibussowitsch return DMPolytopeTypeComposeOrientation(ct, o1, triInv[o2 + 3]); 1108d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUADRILATERAL: 1109d71ae5a4SJacob Faibussowitsch return DMPolytopeTypeComposeOrientation(ct, o1, quadInv[o2 + 4]); 1110d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_SEG_PRISM_TENSOR: 1111d71ae5a4SJacob Faibussowitsch return DMPolytopeTypeComposeOrientation(ct, o1, o2); 1112d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TETRAHEDRON: 1113d71ae5a4SJacob Faibussowitsch return DMPolytopeTypeComposeOrientation(ct, o1, tetInv[o2 + 12]); 1114d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_HEXAHEDRON: 1115d71ae5a4SJacob Faibussowitsch return DMPolytopeTypeComposeOrientation(ct, o1, hexInv[o2 + 24]); 1116d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRI_PRISM: 1117d71ae5a4SJacob Faibussowitsch return DMPolytopeTypeComposeOrientation(ct, o1, tripInv[o2 + 6]); 1118d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRI_PRISM_TENSOR: 1119d71ae5a4SJacob Faibussowitsch return DMPolytopeTypeComposeOrientation(ct, o1, ttriInv[o2 + 6]); 1120d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUAD_PRISM_TENSOR: 1121d71ae5a4SJacob Faibussowitsch return DMPolytopeTypeComposeOrientation(ct, o1, tquadInv[o2 + 8]); 1122d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_PYRAMID: 1123d71ae5a4SJacob Faibussowitsch return DMPolytopeTypeComposeOrientation(ct, o1, pyrInv[o2 + 4]); 1124d71ae5a4SJacob Faibussowitsch default: 1125d71ae5a4SJacob Faibussowitsch return 0; 1126b5a892a1SMatthew G. Knepley } 1127b5a892a1SMatthew G. Knepley } 1128b5a892a1SMatthew G. Knepley 1129b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeMatchOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *); 1130b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeMatchVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *); 1131b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeGetOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *); 1132b5a892a1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeGetVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *); 1133012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPolytopeInCellTest(DMPolytopeType, const PetscReal[], PetscBool *); 1134