1 /* 2 Objects to manage the interactions between the mesh data structures and the algebraic objects 3 */ 4 #ifndef PETSCDM_H 5 #define PETSCDM_H 6 #include <petscmat.h> 7 #include <petscdmtypes.h> 8 #include <petscfetypes.h> 9 #include <petscdstypes.h> 10 #include <petscdmlabel.h> 11 #include <petscdt.h> 12 13 /* SUBMANSEC = DM */ 14 15 PETSC_EXTERN PetscErrorCode DMInitializePackage(void); 16 17 PETSC_EXTERN PetscClassId DM_CLASSID; 18 19 #define DMLOCATEPOINT_POINT_NOT_FOUND -367 20 21 /*J 22 DMType - String with the name of a PETSc `DM` 23 24 Level: beginner 25 26 .seealso: [](ch_dmbase), `DMSetType()`, `DMCreate()`, `DM` 27 J*/ 28 typedef const char *DMType; 29 #define DMDA "da" 30 #define DMCOMPOSITE "composite" 31 #define DMSLICED "sliced" 32 #define DMSHELL "shell" 33 #define DMPLEX "plex" 34 #define DMREDUNDANT "redundant" 35 #define DMPATCH "patch" 36 #define DMMOAB "moab" 37 #define DMNETWORK "network" 38 #define DMFOREST "forest" 39 #define DMP4EST "p4est" 40 #define DMP8EST "p8est" 41 #define DMSWARM "swarm" 42 #define DMPRODUCT "product" 43 #define DMSTAG "stag" 44 45 PETSC_EXTERN const char *const DMBoundaryTypes[]; 46 PETSC_EXTERN const char *const DMBoundaryConditionTypes[]; 47 PETSC_EXTERN const char *const DMBlockingTypes[]; 48 PETSC_EXTERN PetscFunctionList DMList; 49 PETSC_EXTERN DMGeneratorFunctionList DMGenerateList; 50 PETSC_EXTERN PetscErrorCode DMCreate(MPI_Comm, DM *); 51 PETSC_EXTERN PetscErrorCode DMClone(DM, DM *); 52 PETSC_EXTERN PetscErrorCode DMSetType(DM, DMType); 53 PETSC_EXTERN PetscErrorCode DMGetType(DM, DMType *); 54 PETSC_EXTERN PetscErrorCode DMRegister(const char[], PetscErrorCode (*)(DM)); 55 PETSC_EXTERN PetscErrorCode DMRegisterDestroy(void); 56 57 PETSC_EXTERN PetscErrorCode DMView(DM, PetscViewer); 58 PETSC_EXTERN PetscErrorCode DMLoad(DM, PetscViewer); 59 PETSC_EXTERN PetscErrorCode DMDestroy(DM *); 60 PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM, Vec *); 61 PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM, Vec *); 62 PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM, Vec *); 63 PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM, Vec *); 64 PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM, Vec *); 65 PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM, Vec *); 66 PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM); 67 PETSC_EXTERN PetscErrorCode DMClearLocalVectors(DM); 68 PETSC_EXTERN PetscErrorCode DMHasNamedGlobalVector(DM, const char *, PetscBool *); 69 PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM, const char *, Vec *); 70 PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM, const char *, Vec *); 71 PETSC_EXTERN PetscErrorCode DMHasNamedLocalVector(DM, const char *, PetscBool *); 72 PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM, const char *, Vec *); 73 PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM, const char *, Vec *); 74 PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM, ISLocalToGlobalMapping *); 75 PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM, PetscInt *, char ***, IS **); 76 PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM, PetscInt *); 77 PETSC_EXTERN PetscErrorCode DMCreateColoring(DM, ISColoringType, ISColoring *); 78 PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM, Mat *); 79 PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateSkip(DM, PetscBool); 80 PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM, PetscBool); 81 PETSC_EXTERN PetscErrorCode DMSetMatrixStructureOnly(DM, PetscBool); 82 PETSC_EXTERN PetscErrorCode DMSetBlockingType(DM, DMBlockingType); 83 PETSC_EXTERN PetscErrorCode DMGetBlockingType(DM, DMBlockingType *); 84 PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM, DM, Mat *, Vec *); 85 PETSC_EXTERN PetscErrorCode DMCreateRestriction(DM, DM, Mat *); 86 PETSC_EXTERN PetscErrorCode DMCreateInjection(DM, DM, Mat *); 87 PETSC_EXTERN PetscErrorCode DMCreateMassMatrix(DM, DM, Mat *); 88 PETSC_EXTERN PetscErrorCode DMCreateMassMatrixLumped(DM, Vec *); 89 PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM, PetscInt, MPI_Datatype, void *); 90 PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM, PetscInt, MPI_Datatype, void *); 91 PETSC_EXTERN PetscErrorCode DMRefine(DM, MPI_Comm, DM *); 92 PETSC_EXTERN PetscErrorCode DMCoarsen(DM, MPI_Comm, DM *); 93 PETSC_EXTERN PetscErrorCode DMGetCoarseDM(DM, DM *); 94 PETSC_EXTERN PetscErrorCode DMSetCoarseDM(DM, DM); 95 PETSC_EXTERN PetscErrorCode DMGetFineDM(DM, DM *); 96 PETSC_EXTERN PetscErrorCode DMSetFineDM(DM, DM); 97 PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM, PetscInt, DM[]); 98 PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM, PetscInt, DM[]); 99 PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, Vec, Mat, DM, void *), void *); 100 PETSC_EXTERN PetscErrorCode DMCoarsenHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, Vec, Mat, DM, void *), void *); 101 PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, DM, void *), void *); 102 PETSC_EXTERN PetscErrorCode DMRefineHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, Mat, DM, void *), void *); 103 PETSC_EXTERN PetscErrorCode DMRestrict(DM, Mat, Vec, Mat, DM); 104 PETSC_EXTERN PetscErrorCode DMInterpolate(DM, Mat, DM); 105 PETSC_EXTERN PetscErrorCode DMInterpolateSolution(DM, DM, Mat, Vec, Vec); 106 PETSC_EXTERN PetscErrorCode DMExtrude(DM, PetscInt, DM *); 107 PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM); 108 PETSC_EXTERN PetscErrorCode DMViewFromOptions(DM, PetscObject, const char[]); 109 110 PETSC_EXTERN PetscErrorCode DMGenerate(DM, const char[], PetscBool, DM *); 111 PETSC_EXTERN PetscErrorCode DMGenerateRegister(const char[], PetscErrorCode (*)(DM, PetscBool, DM *), PetscErrorCode (*)(DM, PetscReal *, DM *), PetscErrorCode (*)(DM, Vec, DMLabel, DMLabel, DM *), PetscInt); 112 PETSC_EXTERN PetscErrorCode DMGenerateRegisterAll(void); 113 PETSC_EXTERN PetscErrorCode DMGenerateRegisterDestroy(void); 114 PETSC_EXTERN PetscErrorCode DMAdaptLabel(DM, DMLabel, DM *); 115 PETSC_EXTERN PetscErrorCode DMAdaptMetric(DM, Vec, DMLabel, DMLabel, DM *); 116 117 PETSC_EXTERN PetscErrorCode DMSetUp(DM); 118 PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM, DM, Mat, Vec *); 119 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use DMDACreateAggregates() or DMCreateRestriction() (since version 3.12)") PetscErrorCode DMCreateAggregates(DM, DM, Mat *); 120 PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM, PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), void *); 121 PETSC_EXTERN PetscErrorCode DMLocalToGlobalHookAdd(DM, PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), PetscErrorCode (*)(DM, Vec, InsertMode, Vec, void *), void *); 122 PETSC_EXTERN PetscErrorCode DMGlobalToLocal(DM, Vec, InsertMode, Vec); 123 PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM, Vec, InsertMode, Vec); 124 PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM, Vec, InsertMode, Vec); 125 PETSC_EXTERN PetscErrorCode DMLocalToGlobal(DM, Vec, InsertMode, Vec); 126 PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM, Vec, InsertMode, Vec); 127 PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM, Vec, InsertMode, Vec); 128 PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM, Vec, InsertMode, Vec); 129 PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM, Vec, InsertMode, Vec); 130 PETSC_EXTERN PetscErrorCode DMConvert(DM, DMType, DM *); 131 132 /* Topology support */ 133 PETSC_EXTERN PetscErrorCode DMGetDimension(DM, PetscInt *); 134 PETSC_EXTERN PetscErrorCode DMSetDimension(DM, PetscInt); 135 PETSC_EXTERN PetscErrorCode DMGetDimPoints(DM, PetscInt, PetscInt *, PetscInt *); 136 PETSC_EXTERN PetscErrorCode DMGetUseNatural(DM, PetscBool *); 137 PETSC_EXTERN PetscErrorCode DMSetUseNatural(DM, PetscBool); 138 PETSC_EXTERN PetscErrorCode DMGetNeighbors(DM, PetscInt *, const PetscMPIInt **); 139 140 /* Coordinate support */ 141 PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM, DM *); 142 PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM, DM); 143 PETSC_EXTERN PetscErrorCode DMGetCellCoordinateDM(DM, DM *); 144 PETSC_EXTERN PetscErrorCode DMSetCellCoordinateDM(DM, DM); 145 PETSC_EXTERN PetscErrorCode DMGetCoordinateDim(DM, PetscInt *); 146 PETSC_EXTERN PetscErrorCode DMSetCoordinateDim(DM, PetscInt); 147 PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM, PetscSection *); 148 PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM, PetscInt, PetscSection); 149 PETSC_EXTERN PetscErrorCode DMGetCellCoordinateSection(DM, PetscSection *); 150 PETSC_EXTERN PetscErrorCode DMSetCellCoordinateSection(DM, PetscInt, PetscSection); 151 PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM, Vec *); 152 PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM, Vec); 153 PETSC_EXTERN PetscErrorCode DMGetCellCoordinates(DM, Vec *); 154 PETSC_EXTERN PetscErrorCode DMSetCellCoordinates(DM, Vec); 155 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalSetUp(DM); 156 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM, Vec *); 157 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalNoncollective(DM, Vec *); 158 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalTuple(DM, IS, PetscSection *, Vec *); 159 PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM, Vec); 160 PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM); 161 PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocal(DM, Vec *); 162 PETSC_EXTERN PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM, Vec *); 163 PETSC_EXTERN PetscErrorCode DMSetCellCoordinatesLocal(DM, Vec); 164 PETSC_EXTERN PetscErrorCode DMGetCoordinateField(DM, DMField *); 165 PETSC_EXTERN PetscErrorCode DMSetCoordinateField(DM, DMField); 166 PETSC_EXTERN PetscErrorCode DMGetLocalBoundingBox(DM, PetscReal[], PetscReal[]); 167 PETSC_EXTERN PetscErrorCode DMGetBoundingBox(DM, PetscReal[], PetscReal[]); 168 PETSC_EXTERN PetscErrorCode DMProjectCoordinates(DM, PetscFE); 169 PETSC_EXTERN PetscErrorCode DMLocatePoints(DM, Vec, DMPointLocationType, PetscSF *); 170 171 /* Periodicity support */ 172 PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM, const PetscReal *[], const PetscReal *[], const PetscReal *[]); 173 PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM, const PetscReal[], const PetscReal[], const PetscReal[]); 174 PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate(DM, const PetscScalar[], PetscBool, PetscScalar[]); 175 PETSC_EXTERN PetscErrorCode DMLocalizeCoordinates(DM); 176 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalized(DM, PetscBool *); 177 PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalizedLocal(DM, PetscBool *); 178 179 /* block hook interface */ 180 PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, VecScatter, VecScatter, DM, void *), void *); 181 PETSC_EXTERN PetscErrorCode DMSubDomainHookRemove(DM, PetscErrorCode (*)(DM, DM, void *), PetscErrorCode (*)(DM, VecScatter, VecScatter, DM, void *), void *); 182 PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM, VecScatter, VecScatter, DM); 183 184 PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM, const char[]); 185 PETSC_EXTERN PetscErrorCode DMAppendOptionsPrefix(DM, const char[]); 186 PETSC_EXTERN PetscErrorCode DMGetOptionsPrefix(DM, const char *[]); 187 PETSC_EXTERN PetscErrorCode DMSetVecType(DM, VecType); 188 PETSC_EXTERN PetscErrorCode DMGetVecType(DM, VecType *); 189 PETSC_EXTERN PetscErrorCode DMSetMatType(DM, MatType); 190 PETSC_EXTERN PetscErrorCode DMGetMatType(DM, MatType *); 191 PETSC_EXTERN PetscErrorCode DMSetISColoringType(DM, ISColoringType); 192 PETSC_EXTERN PetscErrorCode DMGetISColoringType(DM, ISColoringType *); 193 PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM, void *); 194 PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM, PetscErrorCode (*)(void **)); 195 PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM, void *); 196 PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM, PetscErrorCode (*)(DM, Vec, Vec)); 197 PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM, PetscBool *); 198 PETSC_EXTERN PetscErrorCode DMHasColoring(DM, PetscBool *); 199 PETSC_EXTERN PetscErrorCode DMHasCreateRestriction(DM, PetscBool *); 200 PETSC_EXTERN PetscErrorCode DMHasCreateInjection(DM, PetscBool *); 201 PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM, Vec, Vec); 202 203 PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, const PetscInt[], IS *, DM *); 204 PETSC_EXTERN PetscErrorCode DMCreateSuperDM(DM[], PetscInt, IS **, DM *); 205 PETSC_EXTERN PetscErrorCode DMCreateSectionSubDM(DM, PetscInt, const PetscInt[], IS *, DM *); 206 PETSC_EXTERN PetscErrorCode DMCreateSectionSuperDM(DM[], PetscInt, IS **, DM *); 207 PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM, PetscInt *, char ***, IS **, DM **); 208 PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM, PetscInt *, char ***, IS **, IS **, DM **); 209 PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM, PetscInt, DM *, VecScatter **, VecScatter **, VecScatter **); 210 211 PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM, PetscInt *); 212 PETSC_EXTERN PetscErrorCode DMSetRefineLevel(DM, PetscInt); 213 PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM, PetscInt *); 214 PETSC_EXTERN PetscErrorCode DMSetCoarsenLevel(DM, PetscInt); 215 PETSC_EXTERN PetscErrorCode DMFinalizePackage(void); 216 217 PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM *); 218 PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM); 219 PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM *); 220 PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM); 221 PETSC_EXTERN PetscErrorCode MatFDColoringUseDM(Mat, MatFDColoring); 222 223 typedef struct NLF_DAAD *NLF; 224 225 #define DM_FILE_CLASSID 1211221 226 227 /* FEM support */ 228 PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char[], PetscInt, const PetscScalar[]); 229 PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char[], PetscInt, PetscInt, const PetscScalar[]); 230 PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char[], PetscReal, Vec); 231 232 PETSC_EXTERN PetscErrorCode DMSetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *)); 233 PETSC_EXTERN PetscErrorCode DMGetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *)); 234 PETSC_EXTERN PetscErrorCode DMSetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *)); 235 PETSC_EXTERN PetscErrorCode DMGetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *)); 236 237 PETSC_EXTERN PetscErrorCode DMGetSection(DM, PetscSection *); /* Use DMGetLocalSection() in new code (since v3.12) */ 238 PETSC_EXTERN PetscErrorCode DMSetSection(DM, PetscSection); /* Use DMSetLocalSection() in new code (since v3.12) */ 239 PETSC_EXTERN PetscErrorCode DMGetLocalSection(DM, PetscSection *); 240 PETSC_EXTERN PetscErrorCode DMSetLocalSection(DM, PetscSection); 241 PETSC_EXTERN PetscErrorCode DMGetGlobalSection(DM, PetscSection *); 242 PETSC_EXTERN PetscErrorCode DMSetGlobalSection(DM, PetscSection); 243 static inline PETSC_DEPRECATED_FUNCTION("Use DMGetSection() (since v3.9)") PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *s) 244 { 245 return DMGetSection(dm, s); 246 } 247 static inline PETSC_DEPRECATED_FUNCTION("Use DMSetSection() (since v3.9)") PetscErrorCode DMSetDefaultSection(DM dm, PetscSection s) 248 { 249 return DMSetSection(dm, s); 250 } 251 static inline PETSC_DEPRECATED_FUNCTION("Use DMGetGlobalSection() (since v3.9)") PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *s) 252 { 253 return DMGetGlobalSection(dm, s); 254 } 255 static inline PETSC_DEPRECATED_FUNCTION("Use DMSetGlobalSection() (since v3.9)") PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection s) 256 { 257 return DMSetGlobalSection(dm, s); 258 } 259 260 PETSC_EXTERN PetscErrorCode DMGetSectionSF(DM, PetscSF *); 261 PETSC_EXTERN PetscErrorCode DMSetSectionSF(DM, PetscSF); 262 PETSC_EXTERN PetscErrorCode DMCreateSectionSF(DM, PetscSection, PetscSection); 263 static inline PETSC_DEPRECATED_FUNCTION("Use DMGetSectionSF() (since v3.12)") PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *s) 264 { 265 return DMGetSectionSF(dm, s); 266 } 267 static inline PETSC_DEPRECATED_FUNCTION("Use DMSetSectionSF() (since v3.12)") PetscErrorCode DMSetDefaultSF(DM dm, PetscSF s) 268 { 269 return DMSetSectionSF(dm, s); 270 } 271 static inline PETSC_DEPRECATED_FUNCTION("Use DMCreateSectionSF() (since v3.12)") PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection l, PetscSection g) 272 { 273 return DMCreateSectionSF(dm, l, g); 274 } 275 PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *); 276 PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF); 277 PETSC_EXTERN PetscErrorCode DMGetNaturalSF(DM, PetscSF *); 278 PETSC_EXTERN PetscErrorCode DMSetNaturalSF(DM, PetscSF); 279 280 PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *, Vec *); 281 PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat, Vec); 282 283 PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *); 284 PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *); 285 PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal); 286 PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char *, PetscInt, PetscReal *); 287 288 PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *); 289 PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt); 290 PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, DMLabel *, PetscObject *); 291 PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, DMLabel, PetscObject); 292 PETSC_EXTERN PetscErrorCode DMAddField(DM, DMLabel, PetscObject); 293 PETSC_EXTERN PetscErrorCode DMSetFieldAvoidTensor(DM, PetscInt, PetscBool); 294 PETSC_EXTERN PetscErrorCode DMGetFieldAvoidTensor(DM, PetscInt, PetscBool *); 295 PETSC_EXTERN PetscErrorCode DMClearFields(DM); 296 PETSC_EXTERN PetscErrorCode DMCopyFields(DM, DM); 297 PETSC_EXTERN PetscErrorCode DMGetAdjacency(DM, PetscInt, PetscBool *, PetscBool *); 298 PETSC_EXTERN PetscErrorCode DMSetAdjacency(DM, PetscInt, PetscBool, PetscBool); 299 PETSC_EXTERN PetscErrorCode DMGetBasicAdjacency(DM, PetscBool *, PetscBool *); 300 PETSC_EXTERN PetscErrorCode DMSetBasicAdjacency(DM, PetscBool, PetscBool); 301 302 PETSC_EXTERN PetscErrorCode DMGetNumDS(DM, PetscInt *); 303 PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *); 304 PETSC_EXTERN PetscErrorCode DMGetCellDS(DM, PetscInt, PetscDS *, PetscDS *); 305 PETSC_EXTERN PetscErrorCode DMGetRegionDS(DM, DMLabel, IS *, PetscDS *, PetscDS *); 306 PETSC_EXTERN PetscErrorCode DMSetRegionDS(DM, DMLabel, IS, PetscDS, PetscDS); 307 PETSC_EXTERN PetscErrorCode DMGetRegionNumDS(DM, PetscInt, DMLabel *, IS *, PetscDS *, PetscDS *); 308 PETSC_EXTERN PetscErrorCode DMSetRegionNumDS(DM, PetscInt, DMLabel, IS, PetscDS, PetscDS); 309 PETSC_EXTERN PetscErrorCode DMFindRegionNum(DM, PetscDS, PetscInt *); 310 PETSC_EXTERN PetscErrorCode DMCreateFEDefault(DM, PetscInt, const char[], PetscInt, PetscFE *); 311 PETSC_EXTERN PetscErrorCode DMCreateDS(DM); 312 PETSC_EXTERN PetscErrorCode DMClearDS(DM); 313 PETSC_EXTERN PetscErrorCode DMCopyDS(DM, DM); 314 PETSC_EXTERN PetscErrorCode DMCopyDisc(DM, DM); 315 PETSC_EXTERN PetscErrorCode DMComputeExactSolution(DM, PetscReal, Vec, Vec); 316 PETSC_EXTERN PetscErrorCode DMGetNumAuxiliaryVec(DM, PetscInt *); 317 PETSC_EXTERN PetscErrorCode DMGetAuxiliaryVec(DM, DMLabel, PetscInt, PetscInt, Vec *); 318 PETSC_EXTERN PetscErrorCode DMSetAuxiliaryVec(DM, DMLabel, PetscInt, PetscInt, Vec); 319 PETSC_EXTERN PetscErrorCode DMGetAuxiliaryLabels(DM, DMLabel[], PetscInt[], PetscInt[]); 320 PETSC_EXTERN PetscErrorCode DMCopyAuxiliaryVec(DM, DM); 321 322 /*MC 323 DMInterpolationInfo - Structure for holding information about interpolation on a mesh 324 325 Synopsis: 326 comm - The communicator 327 dim - The spatial dimension of points 328 nInput - The number of input points 329 points - The input point coordinates 330 cells - The cell containing each point 331 n - The number of local points 332 coords - The point coordinates 333 dof - The number of components to interpolate 334 335 Level: intermediate 336 337 .seealso: `DM`, `DMInterpolationCreate()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()` 338 M*/ 339 struct _DMInterpolationInfo { 340 MPI_Comm comm; 341 PetscInt dim; /*1 The spatial dimension of points */ 342 PetscInt nInput; /* The number of input points */ 343 PetscReal *points; /* The input point coordinates */ 344 PetscInt *cells; /* The cell containing each point */ 345 PetscInt n; /* The number of local points */ 346 Vec coords; /* The point coordinates */ 347 PetscInt dof; /* The number of components to interpolate */ 348 }; 349 typedef struct _DMInterpolationInfo *DMInterpolationInfo; 350 351 PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *); 352 PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt); 353 PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *); 354 PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt); 355 PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *); 356 PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]); 357 PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool, PetscBool); 358 PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *); 359 PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *); 360 PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *); 361 PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec); 362 PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *); 363 364 PETSC_EXTERN PetscErrorCode DMCreateLabel(DM, const char[]); 365 PETSC_EXTERN PetscErrorCode DMGetLabelValue(DM, const char[], PetscInt, PetscInt *); 366 PETSC_EXTERN PetscErrorCode DMSetLabelValue(DM, const char[], PetscInt, PetscInt); 367 PETSC_EXTERN PetscErrorCode DMClearLabelValue(DM, const char[], PetscInt, PetscInt); 368 PETSC_EXTERN PetscErrorCode DMGetLabelSize(DM, const char[], PetscInt *); 369 PETSC_EXTERN PetscErrorCode DMGetLabelIdIS(DM, const char[], IS *); 370 PETSC_EXTERN PetscErrorCode DMGetStratumSize(DM, const char[], PetscInt, PetscInt *); 371 PETSC_EXTERN PetscErrorCode DMGetStratumIS(DM, const char[], PetscInt, IS *); 372 PETSC_EXTERN PetscErrorCode DMSetStratumIS(DM, const char[], PetscInt, IS); 373 PETSC_EXTERN PetscErrorCode DMClearLabelStratum(DM, const char[], PetscInt); 374 PETSC_EXTERN PetscErrorCode DMGetLabelOutput(DM, const char[], PetscBool *); 375 PETSC_EXTERN PetscErrorCode DMSetLabelOutput(DM, const char[], PetscBool); 376 PETSC_EXTERN PetscErrorCode DMGetFirstLabeledPoint(DM, DM, DMLabel, PetscInt, const PetscInt *, PetscInt, PetscInt *, PetscDS *); 377 378 /*E 379 DMCopyLabelsMode - Determines how `DMCopyLabels()` behaves when there is a `DMLabel` in the source and destination DMs with the same name 380 381 Values: 382 + `DM_COPY_LABELS_REPLACE` - replace label in destination by label from source 383 . `DM_COPY_LABELS_KEEP` - keep destination label 384 - `DM_COPY_LABELS_FAIL` - throw error 385 386 Level: advanced 387 388 .seealso: `DMLabel`, `DM`, `DMCompareLabels()`, `DMRemoveLabel()` 389 E*/ 390 typedef enum { 391 DM_COPY_LABELS_REPLACE, 392 DM_COPY_LABELS_KEEP, 393 DM_COPY_LABELS_FAIL 394 } DMCopyLabelsMode; 395 PETSC_EXTERN const char *const DMCopyLabelsModes[]; 396 397 PETSC_EXTERN PetscErrorCode DMGetNumLabels(DM, PetscInt *); 398 PETSC_EXTERN PetscErrorCode DMGetLabelName(DM, PetscInt, const char **); 399 PETSC_EXTERN PetscErrorCode DMHasLabel(DM, const char[], PetscBool *); 400 PETSC_EXTERN PetscErrorCode DMGetLabel(DM, const char *, DMLabel *); 401 PETSC_EXTERN PetscErrorCode DMSetLabel(DM, DMLabel); 402 PETSC_EXTERN PetscErrorCode DMGetLabelByNum(DM, PetscInt, DMLabel *); 403 PETSC_EXTERN PetscErrorCode DMAddLabel(DM, DMLabel); 404 PETSC_EXTERN PetscErrorCode DMRemoveLabel(DM, const char[], DMLabel *); 405 PETSC_EXTERN PetscErrorCode DMRemoveLabelBySelf(DM, DMLabel *, PetscBool); 406 PETSC_EXTERN PetscErrorCode DMCopyLabels(DM, DM, PetscCopyMode, PetscBool, DMCopyLabelsMode emode); 407 PETSC_EXTERN PetscErrorCode DMCompareLabels(DM, DM, PetscBool *, char **); 408 409 PETSC_EXTERN PetscErrorCode DMAddBoundary(DM, DMBoundaryConditionType, const char[], DMLabel, PetscInt, const PetscInt[], PetscInt, PetscInt, const PetscInt[], void (*)(void), void (*)(void), void *, PetscInt *); 410 PETSC_EXTERN PetscErrorCode DMIsBoundaryPoint(DM, PetscInt, PetscBool *); 411 412 PETSC_EXTERN PetscErrorCode DMProjectFunction(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec); 413 PETSC_EXTERN PetscErrorCode DMProjectFunctionLocal(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec); 414 PETSC_EXTERN PetscErrorCode DMProjectFunctionLabel(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec); 415 PETSC_EXTERN PetscErrorCode DMProjectFunctionLabelLocal(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec); 416 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); 417 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); 418 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); 419 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); 420 PETSC_EXTERN PetscErrorCode DMComputeL2Diff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *); 421 PETSC_EXTERN PetscErrorCode DMComputeL2GradientDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal[], PetscReal *); 422 PETSC_EXTERN PetscErrorCode DMComputeL2FieldDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **, Vec, PetscReal *); 423 PETSC_EXTERN PetscErrorCode DMComputeError(DM, Vec, PetscReal[], Vec *); 424 PETSC_EXTERN PetscErrorCode DMHasBasisTransform(DM, PetscBool *); 425 PETSC_EXTERN PetscErrorCode DMCopyTransform(DM, DM); 426 427 PETSC_EXTERN PetscErrorCode DMGetCompatibility(DM, DM, PetscBool *, PetscBool *); 428 429 PETSC_EXTERN PetscErrorCode DMMonitorSet(DM, PetscErrorCode (*)(DM, void *), void *, PetscErrorCode (*)(void **)); 430 PETSC_EXTERN PetscErrorCode DMMonitorCancel(DM); 431 PETSC_EXTERN PetscErrorCode DMMonitorSetFromOptions(DM, const char[], const char[], const char[], PetscErrorCode (*)(DM, void *), PetscErrorCode (*)(DM, PetscViewerAndFormat *), PetscBool *); 432 PETSC_EXTERN PetscErrorCode DMMonitor(DM); 433 434 static inline PetscInt DMPolytopeTypeGetDim(DMPolytopeType ct) 435 { 436 switch (ct) { 437 case DM_POLYTOPE_POINT: 438 return 0; 439 case DM_POLYTOPE_SEGMENT: 440 case DM_POLYTOPE_POINT_PRISM_TENSOR: 441 return 1; 442 case DM_POLYTOPE_TRIANGLE: 443 case DM_POLYTOPE_QUADRILATERAL: 444 case DM_POLYTOPE_SEG_PRISM_TENSOR: 445 return 2; 446 case DM_POLYTOPE_TETRAHEDRON: 447 case DM_POLYTOPE_HEXAHEDRON: 448 case DM_POLYTOPE_TRI_PRISM: 449 case DM_POLYTOPE_TRI_PRISM_TENSOR: 450 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 451 case DM_POLYTOPE_PYRAMID: 452 return 3; 453 default: 454 return -1; 455 } 456 } 457 458 static inline PetscInt DMPolytopeTypeGetConeSize(DMPolytopeType ct) 459 { 460 switch (ct) { 461 case DM_POLYTOPE_POINT: 462 return 0; 463 case DM_POLYTOPE_SEGMENT: 464 return 2; 465 case DM_POLYTOPE_POINT_PRISM_TENSOR: 466 return 2; 467 case DM_POLYTOPE_TRIANGLE: 468 return 3; 469 case DM_POLYTOPE_QUADRILATERAL: 470 return 4; 471 case DM_POLYTOPE_SEG_PRISM_TENSOR: 472 return 4; 473 case DM_POLYTOPE_TETRAHEDRON: 474 return 4; 475 case DM_POLYTOPE_HEXAHEDRON: 476 return 6; 477 case DM_POLYTOPE_TRI_PRISM: 478 return 5; 479 case DM_POLYTOPE_TRI_PRISM_TENSOR: 480 return 5; 481 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 482 return 6; 483 case DM_POLYTOPE_PYRAMID: 484 return 5; 485 default: 486 return -1; 487 } 488 } 489 490 static inline PetscInt DMPolytopeTypeGetNumVertices(DMPolytopeType ct) 491 { 492 switch (ct) { 493 case DM_POLYTOPE_POINT: 494 return 1; 495 case DM_POLYTOPE_SEGMENT: 496 return 2; 497 case DM_POLYTOPE_POINT_PRISM_TENSOR: 498 return 2; 499 case DM_POLYTOPE_TRIANGLE: 500 return 3; 501 case DM_POLYTOPE_QUADRILATERAL: 502 return 4; 503 case DM_POLYTOPE_SEG_PRISM_TENSOR: 504 return 4; 505 case DM_POLYTOPE_TETRAHEDRON: 506 return 4; 507 case DM_POLYTOPE_HEXAHEDRON: 508 return 8; 509 case DM_POLYTOPE_TRI_PRISM: 510 return 6; 511 case DM_POLYTOPE_TRI_PRISM_TENSOR: 512 return 6; 513 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 514 return 8; 515 case DM_POLYTOPE_PYRAMID: 516 return 5; 517 default: 518 return -1; 519 } 520 } 521 522 static inline DMPolytopeType DMPolytopeTypeSimpleShape(PetscInt dim, PetscBool simplex) 523 { 524 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))); 525 } 526 527 static inline PetscInt DMPolytopeTypeGetNumArrangments(DMPolytopeType ct) 528 { 529 switch (ct) { 530 case DM_POLYTOPE_POINT: 531 return 1; 532 case DM_POLYTOPE_SEGMENT: 533 return 2; 534 case DM_POLYTOPE_POINT_PRISM_TENSOR: 535 return 2; 536 case DM_POLYTOPE_TRIANGLE: 537 return 6; 538 case DM_POLYTOPE_QUADRILATERAL: 539 return 8; 540 case DM_POLYTOPE_SEG_PRISM_TENSOR: 541 return 4; 542 case DM_POLYTOPE_TETRAHEDRON: 543 return 24; 544 case DM_POLYTOPE_HEXAHEDRON: 545 return 48; 546 case DM_POLYTOPE_TRI_PRISM: 547 return 12; 548 case DM_POLYTOPE_TRI_PRISM_TENSOR: 549 return 12; 550 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 551 return 16; 552 case DM_POLYTOPE_PYRAMID: 553 return 8; 554 default: 555 return -1; 556 } 557 } 558 559 /* An arrangement is a face order combined with an orientation for each face */ 560 static inline const PetscInt *DMPolytopeTypeGetArrangment(DMPolytopeType ct, PetscInt o) 561 { 562 static const PetscInt pntArr[1 * 2] = {0, 0}; 563 /* a: swap */ 564 static const PetscInt segArr[2 * 2 * 2] = {1, 0, 0, 0, /* -1: a */ 565 0, 0, 1, 0, 566 /* 0: e */}; 567 /* a: swap first two 568 b: swap last two */ 569 static const PetscInt triArr[6 * 3 * 2] = {0, -1, 2, -1, 1, -1, /* -3: b */ 570 2, -1, 1, -1, 0, -1, /* -2: aba */ 571 1, -1, 0, -1, 2, -1, /* -1: a */ 572 0, 0, 1, 0, 2, 0, /* 0: identity */ 573 1, 0, 2, 0, 0, 0, /* 1: ba */ 574 2, 0, 0, 0, 1, 0, 575 /* 2: ab */}; 576 /* a: forward cyclic permutation 577 b: swap first and last pairs */ 578 static const PetscInt quadArr[8 * 4 * 2] = {1, -1, 0, -1, 3, -1, 2, -1, /* -4: b */ 579 0, -1, 3, -1, 2, -1, 1, -1, /* -3: b a^3 = a b */ 580 3, -1, 2, -1, 1, -1, 0, -1, /* -2: b a^2 = a^2 b */ 581 2, -1, 1, -1, 0, -1, 3, -1, /* -1: b a = a^3 b */ 582 0, 0, 1, 0, 2, 0, 3, 0, /* 0: identity */ 583 1, 0, 2, 0, 3, 0, 0, 0, /* 1: a */ 584 2, 0, 3, 0, 0, 0, 1, 0, /* 2: a^2 */ 585 3, 0, 0, 0, 1, 0, 2, 0, 586 /* 3: a^3 */}; 587 /* r: rotate 180 588 b: swap top and bottom segments */ 589 static const PetscInt tsegArr[4 * 4 * 2] = {1, -1, 0, -1, 3, -1, 2, -1, /* -2: r b */ 590 0, -1, 1, -1, 3, 0, 2, 0, /* -1: r */ 591 0, 0, 1, 0, 2, 0, 3, 0, /* 0: identity */ 592 1, 0, 0, 0, 2, -1, 3, -1, 593 /* 1: b */}; 594 /* https://en.wikiversity.org/wiki/Symmetric_group_S4 */ 595 static const PetscInt tetArr[24 * 4 * 2] = {3, -2, 2, -3, 0, -1, 1, -1, /* -12: (1324) p22 */ 596 3, -1, 1, -3, 2, -1, 0, -1, /* -11: (14) p21 */ 597 3, -3, 0, -3, 1, -1, 2, -1, /* -10: (1234) p18 */ 598 2, -1, 3, -1, 1, -3, 0, -2, /* -9: (1423) p17 */ 599 2, -3, 0, -1, 3, -2, 1, -3, /* -8: (1342) p13 */ 600 2, -2, 1, -2, 0, -2, 3, -2, /* -7: (24) p14 */ 601 1, -2, 0, -2, 2, -2, 3, -1, /* -6: (34) p6 */ 602 1, -1, 3, -3, 0, -3, 2, -2, /* -5: (1243) p10 */ 603 1, -3, 2, -1, 3, -1, 0, -3, /* -4: (1432) p9 */ 604 0, -3, 1, -1, 3, -3, 2, -3, /* -3: (12) p1 */ 605 0, -2, 2, -2, 1, -2, 3, -3, /* -2: (23) p2 */ 606 0, -1, 3, -2, 2, -3, 1, -2, /* -1: (13) p5 */ 607 0, 0, 1, 0, 2, 0, 3, 0, /* 0: () p0 */ 608 0, 1, 3, 1, 1, 2, 2, 0, /* 1: (123) p4 */ 609 0, 2, 2, 1, 3, 0, 1, 2, /* 2: (132) p3 */ 610 1, 2, 0, 1, 3, 1, 2, 2, /* 3: (12)(34) p7 */ 611 1, 0, 2, 0, 0, 0, 3, 1, /* 4: (243) p8 */ 612 1, 1, 3, 2, 2, 2, 0, 0, /* 5: (143) p11 */ 613 2, 1, 3, 0, 0, 2, 1, 0, /* 6: (13)(24) p16 */ 614 2, 2, 1, 1, 3, 2, 0, 2, /* 7: (142) p15 */ 615 2, 0, 0, 0, 1, 0, 3, 2, /* 8: (234) p12 */ 616 3, 2, 2, 2, 1, 1, 0, 1, /* 9: (14)(23) p23 */ 617 3, 0, 0, 2, 2, 1, 1, 1, /* 10: (134) p19 */ 618 3, 1, 1, 2, 0, 1, 2, 1 /* 11: (124) p20 */}; 619 /* Each rotation determines a permutation of the four diagonals, and this defines the isomorphism with S_4 */ 620 static const PetscInt hexArr[48 * 6 * 2] = { 621 2, -3, 3, -2, 4, -2, 5, -3, 1, -3, 0, -1, /* -24: reflect bottom and use -3 on top */ 622 4, -2, 5, -2, 0, -1, 1, -4, 3, -2, 2, -3, /* -23: reflect bottom and use -3 on top */ 623 5, -3, 4, -1, 1, -2, 0, -3, 3, -4, 2, -1, /* -22: reflect bottom and use -3 on top */ 624 3, -1, 2, -4, 4, -4, 5, -1, 0, -4, 1, -4, /* -21: reflect bottom and use -3 on top */ 625 3, -3, 2, -2, 5, -1, 4, -4, 1, -1, 0, -3, /* -20: reflect bottom and use -3 on top */ 626 4, -4, 5, -4, 1, -4, 0, -1, 2, -4, 3, -1, /* -19: reflect bottom and use -3 on top */ 627 2, -1, 3, -4, 5, -3, 4, -2, 0, -2, 1, -2, /* -18: reflect bottom and use -3 on top */ 628 5, -1, 4, -3, 0, -3, 1, -2, 2, -2, 3, -3, /* -17: reflect bottom and use -3 on top */ 629 4, -3, 5, -1, 3, -2, 2, -4, 1, -4, 0, -4, /* -16: reflect bottom and use -3 on top */ 630 5, -4, 4, -4, 3, -4, 2, -2, 0, -3, 1, -1, /* -15: reflect bottom and use -3 on top */ 631 3, -4, 2, -1, 1, -1, 0, -4, 4, -4, 5, -4, /* -14: reflect bottom and use -3 on top */ 632 2, -2, 3, -3, 0, -2, 1, -3, 4, -2, 5, -2, /* -13: reflect bottom and use -3 on top */ 633 1, -3, 0, -1, 4, -1, 5, -4, 3, -1, 2, -4, /* -12: reflect bottom and use -3 on top */ 634 1, -1, 0, -3, 5, -4, 4, -1, 2, -1, 3, -4, /* -11: reflect bottom and use -3 on top */ 635 5, -2, 4, -2, 2, -2, 3, -4, 1, -2, 0, -2, /* -10: reflect bottom and use -3 on top */ 636 1, -2, 0, -2, 2, -1, 3, -1, 4, -1, 5, -3, /* -9: reflect bottom and use -3 on top */ 637 4, -1, 5, -3, 2, -4, 3, -2, 0, -1, 1, -3, /* -8: reflect bottom and use -3 on top */ 638 3, -2, 2, -3, 0, -4, 1, -1, 5, -1, 4, -3, /* -7: reflect bottom and use -3 on top */ 639 1, -4, 0, -4, 3, -1, 2, -1, 5, -4, 4, -4, /* -6: reflect bottom and use -3 on top */ 640 2, -4, 3, -1, 1, -3, 0, -2, 5, -3, 4, -1, /* -5: reflect bottom and use -3 on top */ 641 0, -4, 1, -4, 4, -3, 5, -2, 2, -3, 3, -2, /* -4: reflect bottom and use -3 on top */ 642 0, -3, 1, -1, 3, -3, 2, -3, 4, -3, 5, -1, /* -3: reflect bottom and use -3 on top */ 643 0, -2, 1, -2, 5, -2, 4, -3, 3, -3, 2, -2, /* -2: reflect bottom and use -3 on top */ 644 0, -1, 1, -3, 2, -3, 3, -3, 5, -2, 4, -2, /* -1: reflect bottom and use -3 on top */ 645 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, 5, 0, /* 0: identity */ 646 0, 1, 1, 3, 5, 3, 4, 0, 2, 0, 3, 1, /* 1: 90 rotation about z */ 647 0, 2, 1, 2, 3, 0, 2, 0, 5, 3, 4, 1, /* 2: 180 rotation about z */ 648 0, 3, 1, 1, 4, 0, 5, 3, 3, 0, 2, 1, /* 3: 270 rotation about z */ 649 2, 3, 3, 2, 1, 0, 0, 3, 4, 3, 5, 1, /* 4: 90 rotation about x */ 650 1, 3, 0, 1, 3, 2, 2, 2, 4, 2, 5, 2, /* 5: 180 rotation about x */ 651 3, 1, 2, 0, 0, 1, 1, 2, 4, 1, 5, 3, /* 6: 270 rotation about x */ 652 4, 0, 5, 0, 2, 1, 3, 3, 1, 1, 0, 3, /* 7: 90 rotation about y */ 653 1, 1, 0, 3, 2, 2, 3, 2, 5, 1, 4, 3, /* 8: 180 rotation about y */ 654 5, 1, 4, 3, 2, 3, 3, 1, 0, 0, 1, 0, /* 9: 270 rotation about y */ 655 1, 0, 0, 0, 5, 1, 4, 2, 3, 2, 2, 3, /* 10: 180 rotation about x+y */ 656 1, 2, 0, 2, 4, 2, 5, 1, 2, 2, 3, 3, /* 11: 180 rotation about x-y */ 657 2, 1, 3, 0, 0, 3, 1, 0, 5, 0, 4, 0, /* 12: 180 rotation about y+z */ 658 3, 3, 2, 2, 1, 2, 0, 1, 5, 2, 4, 2, /* 13: 180 rotation about y-z */ 659 5, 3, 4, 1, 3, 1, 2, 3, 1, 3, 0, 1, /* 14: 180 rotation about z+x */ 660 4, 2, 5, 2, 3, 3, 2, 1, 0, 2, 1, 2, /* 15: 180 rotation about z-x */ 661 5, 0, 4, 0, 0, 0, 1, 3, 3, 1, 2, 0, /* 16: 120 rotation about x+y+z (v0v6) */ 662 2, 0, 3, 1, 5, 0, 4, 3, 1, 0, 0, 0, /* 17: 240 rotation about x+y+z (v0v6) */ 663 4, 3, 5, 1, 1, 1, 0, 2, 3, 3, 2, 2, /* 18: 120 rotation about x+y-z (v4v2) */ 664 3, 2, 2, 3, 5, 2, 4, 1, 0, 1, 1, 3, /* 19: 240 rotation about x+y-z (v4v2) */ 665 3, 0, 2, 1, 4, 1, 5, 2, 1, 2, 0, 2, /* 20: 120 rotation about x-y+z (v1v5) */ 666 5, 2, 4, 2, 1, 3, 0, 0, 2, 3, 3, 2, /* 21: 240 rotation about x-y+z (v1v5) */ 667 4, 1, 5, 3, 0, 2, 1, 1, 2, 1, 3, 0, /* 22: 120 rotation about x-y-z (v7v3) */ 668 2, 2, 3, 3, 4, 3, 5, 0, 0, 3, 1, 1, /* 23: 240 rotation about x-y-z (v7v3) */ 669 }; 670 static const PetscInt tripArr[12 * 5 * 2] = { 671 1, -3, 0, -1, 3, -1, 4, -1, 2, -1, /* -6: reflect bottom and top */ 672 1, -1, 0, -3, 4, -1, 2, -1, 3, -1, /* -5: reflect bottom and top */ 673 1, -2, 0, -2, 2, -1, 3, -1, 4, -1, /* -4: reflect bottom and top */ 674 0, -3, 1, -1, 3, -3, 2, -3, 4, -3, /* -3: reflect bottom and top */ 675 0, -2, 1, -2, 4, -3, 3, -3, 2, -3, /* -2: reflect bottom and top */ 676 0, -1, 1, -3, 2, -3, 4, -3, 3, -3, /* -1: reflect bottom and top */ 677 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, /* 0: identity */ 678 0, 1, 1, 2, 4, 0, 2, 0, 3, 0, /* 1: 120 rotation about z */ 679 0, 2, 1, 1, 3, 0, 4, 0, 2, 0, /* 2: 240 rotation about z */ 680 1, 1, 0, 2, 2, 2, 4, 2, 3, 2, /* 3: 180 rotation about y of 0 */ 681 1, 0, 0, 0, 4, 2, 3, 2, 2, 2, /* 4: 180 rotation about y of 1 */ 682 1, 2, 0, 1, 3, 2, 2, 2, 4, 2, /* 5: 180 rotation about y of 2 */ 683 }; 684 /* a: rotate 120 about z 685 b: swap top and bottom segments 686 r: reflect */ 687 static const PetscInt ttriArr[12 * 5 * 2] = { 688 1, -3, 0, -3, 2, -2, 4, -2, 3, -2, /* -6: r b a^2 */ 689 1, -2, 0, -2, 4, -2, 3, -2, 2, -2, /* -5: r b a */ 690 1, -1, 0, -1, 3, -2, 2, -2, 4, -2, /* -4: r b */ 691 0, -3, 1, -3, 2, -1, 4, -1, 3, -1, /* -3: r a^2 */ 692 0, -2, 1, -2, 4, -1, 3, -1, 2, -1, /* -2: r a */ 693 0, -1, 1, -1, 3, -1, 2, -1, 4, -1, /* -1: r */ 694 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, /* 0: identity */ 695 0, 1, 1, 1, 3, 0, 4, 0, 2, 0, /* 1: a */ 696 0, 2, 1, 2, 4, 0, 2, 0, 3, 0, /* 2: a^2 */ 697 1, 0, 0, 0, 2, 1, 3, 1, 4, 1, /* 3: b */ 698 1, 1, 0, 1, 3, 1, 4, 1, 2, 1, /* 4: b a */ 699 1, 2, 0, 2, 4, 1, 2, 1, 3, 1, /* 5: b a^2 */ 700 }; 701 /* a: rotate 90 about z 702 b: swap top and bottom segments 703 r: reflect */ 704 static const PetscInt tquadArr[16 * 6 * 2] = { 705 1, -4, 0, -4, 3, -2, 2, -2, 5, -2, 4, -2, /* -8: r b a^3 */ 706 1, -3, 0, -3, 2, -2, 5, -2, 4, -2, 3, -2, /* -7: r b a^2 */ 707 1, -2, 0, -2, 5, -2, 4, -2, 3, -2, 2, -2, /* -6: r b a */ 708 1, -1, 0, -1, 4, -2, 3, -2, 2, -2, 5, -2, /* -5: r b */ 709 0, -4, 1, -4, 3, -1, 2, -1, 5, -1, 4, -1, /* -4: r a^3 */ 710 0, -3, 1, -3, 2, -1, 5, -1, 4, -1, 3, -1, /* -3: r a^2 */ 711 0, -2, 1, -2, 5, -1, 4, -1, 3, -1, 2, -1, /* -2: r a */ 712 0, -1, 1, -1, 4, -1, 3, -1, 2, -1, 5, -1, /* -1: r */ 713 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, 5, 0, /* 0: identity */ 714 0, 1, 1, 1, 3, 0, 4, 0, 5, 0, 2, 0, /* 1: a */ 715 0, 2, 1, 2, 4, 0, 5, 0, 2, 0, 3, 0, /* 2: a^2 */ 716 0, 3, 1, 3, 5, 0, 2, 0, 3, 0, 4, 0, /* 3: a^3 */ 717 1, 0, 0, 0, 2, 1, 3, 1, 4, 1, 5, 1, /* 4: b */ 718 1, 1, 0, 1, 3, 1, 4, 1, 5, 1, 2, 1, /* 5: b a */ 719 1, 2, 0, 2, 4, 1, 5, 1, 2, 1, 3, 1, /* 6: b a^2 */ 720 1, 3, 0, 3, 5, 1, 2, 1, 3, 1, 4, 1, /* 7: b a^3 */ 721 }; 722 static const PetscInt pyrArr[8 * 5 * 2] = { 723 0, -4, 2, -3, 1, -3, 4, -3, 3, -3, /* -4: Reflect bottom face */ 724 0, -3, 3, -3, 2, -3, 1, -3, 4, -3, /* -3: Reflect bottom face */ 725 0, -2, 4, -3, 3, -3, 2, -3, 1, -3, /* -2: Reflect bottom face */ 726 0, -1, 1, -3, 4, -3, 3, -3, 2, -3, /* -1: Reflect bottom face */ 727 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, /* 0: identity */ 728 0, 1, 4, 0, 1, 0, 2, 0, 3, 0, /* 1: 90 rotation about z */ 729 0, 2, 3, 0, 4, 0, 1, 0, 2, 0, /* 2: 180 rotation about z */ 730 0, 3, 2, 0, 3, 0, 4, 0, 1, 0, /* 3: 270 rotation about z */ 731 }; 732 switch (ct) { 733 case DM_POLYTOPE_POINT: 734 return pntArr; 735 case DM_POLYTOPE_SEGMENT: 736 return &segArr[(o + 1) * 2 * 2]; 737 case DM_POLYTOPE_POINT_PRISM_TENSOR: 738 return &segArr[(o + 1) * 2 * 2]; 739 case DM_POLYTOPE_TRIANGLE: 740 return &triArr[(o + 3) * 3 * 2]; 741 case DM_POLYTOPE_QUADRILATERAL: 742 return &quadArr[(o + 4) * 4 * 2]; 743 case DM_POLYTOPE_SEG_PRISM_TENSOR: 744 return &tsegArr[(o + 2) * 4 * 2]; 745 case DM_POLYTOPE_TETRAHEDRON: 746 return &tetArr[(o + 12) * 4 * 2]; 747 case DM_POLYTOPE_HEXAHEDRON: 748 return &hexArr[(o + 24) * 6 * 2]; 749 case DM_POLYTOPE_TRI_PRISM: 750 return &tripArr[(o + 6) * 5 * 2]; 751 case DM_POLYTOPE_TRI_PRISM_TENSOR: 752 return &ttriArr[(o + 6) * 5 * 2]; 753 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 754 return &tquadArr[(o + 8) * 6 * 2]; 755 case DM_POLYTOPE_PYRAMID: 756 return &pyrArr[(o + 4) * 5 * 2]; 757 default: 758 return PETSC_NULLPTR; 759 } 760 } 761 762 /* A vertex arrangement is a vertex order */ 763 static inline const PetscInt *DMPolytopeTypeGetVertexArrangment(DMPolytopeType ct, PetscInt o) 764 { 765 static const PetscInt pntVerts[1] = {0}; 766 static const PetscInt segVerts[2 * 2] = {1, 0, 0, 1}; 767 static const PetscInt triVerts[6 * 3] = {1, 0, 2, 0, 2, 1, 2, 1, 0, 0, 1, 2, 1, 2, 0, 2, 0, 1}; 768 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}; 769 static const PetscInt tsegVerts[4 * 4] = {3, 2, 1, 0, 1, 0, 3, 2, 0, 1, 2, 3, 2, 3, 0, 1}; 770 static const PetscInt tetVerts[24 * 4] = {2, 3, 1, 0, /* -12: (1324) p22 */ 771 3, 1, 2, 0, /* -11: (14) p21 */ 772 1, 2, 3, 0, /* -10: (1234) p18 */ 773 3, 2, 0, 1, /* -9: (1423) p17 */ 774 2, 0, 3, 1, /* -8: (1342) p13 */ 775 0, 3, 2, 1, /* -7: (24) p14 */ 776 0, 1, 3, 2, /* -6: (34) p6 */ 777 1, 3, 0, 2, /* -5: (1243) p10 */ 778 3, 0, 1, 2, /* -4: (1432 p9 */ 779 1, 0, 2, 3, /* -3: (12) p1 */ 780 0, 2, 1, 3, /* -2: (23) p2 */ 781 2, 1, 0, 3, /* -1: (13) p5 */ 782 0, 1, 2, 3, /* 0: () p0 */ 783 1, 2, 0, 3, /* 1: (123) p4 */ 784 2, 0, 1, 3, /* 2: (132) p3 */ 785 1, 0, 3, 2, /* 3: (12)(34) p7 */ 786 0, 3, 1, 2, /* 4: (243) p8 */ 787 3, 1, 0, 2, /* 5: (143) p11 */ 788 2, 3, 0, 1, /* 6: (13)(24) p16 */ 789 3, 0, 2, 1, /* 7: (142) p15 */ 790 0, 2, 3, 1, /* 8: (234) p12 */ 791 3, 2, 1, 0, /* 9: (14)(23) p23 */ 792 2, 1, 3, 0, /* 10: (134) p19 */ 793 1, 3, 2, 0 /* 11: (124) p20 */}; 794 static const PetscInt hexVerts[48 * 8] = { 795 3, 0, 4, 5, 2, 6, 7, 1, /* -24: reflected 23 */ 796 3, 5, 6, 2, 0, 1, 7, 4, /* -23: reflected 22 */ 797 4, 0, 1, 7, 5, 6, 2, 3, /* -22: reflected 21 */ 798 6, 7, 1, 2, 5, 3, 0, 4, /* -21: reflected 20 */ 799 1, 2, 6, 7, 0, 4, 5, 3, /* -20: reflected 19 */ 800 6, 2, 3, 5, 7, 4, 0, 1, /* -19: reflected 18 */ 801 4, 5, 3, 0, 7, 1, 2, 6, /* -18: reflected 17 */ 802 1, 7, 4, 0, 2, 3, 5, 6, /* -17: reflected 16 */ 803 2, 3, 5, 6, 1, 7, 4, 0, /* -16: reflected 15 */ 804 7, 4, 0, 1, 6, 2, 3, 5, /* -15: reflected 14 */ 805 7, 1, 2, 6, 4, 5, 3, 0, /* -14: reflected 13 */ 806 0, 4, 5, 3, 1, 2, 6, 7, /* -13: reflected 12 */ 807 5, 4, 7, 6, 3, 2, 1, 0, /* -12: reflected 11 */ 808 7, 6, 5, 4, 1, 0, 3, 2, /* -11: reflected 10 */ 809 0, 1, 7, 4, 3, 5, 6, 2, /* -10: reflected 9 */ 810 4, 7, 6, 5, 0, 3, 2, 1, /* -9: reflected 8 */ 811 5, 6, 2, 3, 4, 0, 1, 7, /* -8: reflected 7 */ 812 2, 6, 7, 1, 3, 0, 4, 5, /* -7: reflected 6 */ 813 6, 5, 4, 7, 2, 1, 0, 3, /* -6: reflected 5 */ 814 5, 3, 0, 4, 6, 7, 1, 2, /* -5: reflected 4 */ 815 2, 1, 0, 3, 6, 5, 4, 7, /* -4: reflected 3 */ 816 1, 0, 3, 2, 7, 6, 5, 4, /* -3: reflected 2 */ 817 0, 3, 2, 1, 4, 7, 6, 5, /* -2: reflected 1 */ 818 3, 2, 1, 0, 5, 4, 7, 6, /* -1: reflected 0 */ 819 0, 1, 2, 3, 4, 5, 6, 7, /* 0: identity */ 820 1, 2, 3, 0, 7, 4, 5, 6, /* 1: 90 rotation about z */ 821 2, 3, 0, 1, 6, 7, 4, 5, /* 2: 180 rotation about z */ 822 3, 0, 1, 2, 5, 6, 7, 4, /* 3: 270 rotation about z */ 823 4, 0, 3, 5, 7, 6, 2, 1, /* 4: 90 rotation about x */ 824 7, 4, 5, 6, 1, 2, 3, 0, /* 5: 180 rotation about x */ 825 1, 7, 6, 2, 0, 3, 5, 4, /* 6: 270 rotation about x */ 826 3, 2, 6, 5, 0, 4, 7, 1, /* 7: 90 rotation about y */ 827 5, 6, 7, 4, 3, 0, 1, 2, /* 8: 180 rotation about y */ 828 4, 7, 1, 0, 5, 3, 2, 6, /* 9: 270 rotation about y */ 829 4, 5, 6, 7, 0, 1, 2, 3, /* 10: 180 rotation about x+y */ 830 6, 7, 4, 5, 2, 3, 0, 1, /* 11: 180 rotation about x-y */ 831 3, 5, 4, 0, 2, 1, 7, 6, /* 12: 180 rotation about y+z */ 832 6, 2, 1, 7, 5, 4, 0, 3, /* 13: 180 rotation about y-z */ 833 1, 0, 4, 7, 2, 6, 5, 3, /* 14: 180 rotation about z+x */ 834 6, 5, 3, 2, 7, 1, 0, 4, /* 15: 180 rotation about z-x */ 835 0, 4, 7, 1, 3, 2, 6, 5, /* 16: 120 rotation about x+y+z (v0v6) */ 836 0, 3, 5, 4, 1, 7, 6, 2, /* 17: 240 rotation about x+y+z (v0v6) */ 837 5, 3, 2, 6, 4, 7, 1, 0, /* 18: 120 rotation about x+y-z (v4v2) */ 838 7, 6, 2, 1, 4, 0, 3, 5, /* 19: 240 rotation about x+y-z (v4v2) */ 839 2, 1, 7, 6, 3, 5, 4, 0, /* 20: 120 rotation about x-y+z (v1v5) */ 840 7, 1, 0, 4, 6, 5, 3, 2, /* 21: 240 rotation about x-y+z (v1v5) */ 841 2, 6, 5, 3, 1, 0, 4, 7, /* 22: 120 rotation about x-y-z (v7v3) */ 842 5, 4, 0, 3, 6, 2, 1, 7, /* 23: 240 rotation about x-y-z (v7v3) */ 843 }; 844 static const PetscInt tripVerts[12 * 6] = { 845 4, 3, 5, 2, 1, 0, /* -6: reflect bottom and top */ 846 5, 4, 3, 1, 0, 2, /* -5: reflect bottom and top */ 847 3, 5, 4, 0, 2, 1, /* -4: reflect bottom and top */ 848 1, 0, 2, 5, 4, 3, /* -3: reflect bottom and top */ 849 0, 2, 1, 3, 5, 4, /* -2: reflect bottom and top */ 850 2, 1, 0, 4, 3, 5, /* -1: reflect bottom and top */ 851 0, 1, 2, 3, 4, 5, /* 0: identity */ 852 1, 2, 0, 5, 3, 4, /* 1: 120 rotation about z */ 853 2, 0, 1, 4, 5, 3, /* 2: 240 rotation about z */ 854 4, 5, 3, 2, 0, 1, /* 3: 180 rotation about y of 0 */ 855 3, 4, 5, 0, 1, 2, /* 4: 180 rotation about y of 1 */ 856 5, 3, 4, 1, 2, 0, /* 5: 180 rotation about y of 2 */ 857 }; 858 static const PetscInt ttriVerts[12 * 6] = { 859 4, 3, 5, 1, 0, 2, /* -6: r b a^2 */ 860 3, 5, 4, 0, 2, 1, /* -5: r b a */ 861 5, 4, 3, 2, 1, 0, /* -4: r b */ 862 1, 0, 2, 4, 3, 5, /* -3: r a^2 */ 863 0, 2, 1, 3, 5, 4, /* -2: r a */ 864 2, 1, 0, 5, 4, 3, /* -1: r */ 865 0, 1, 2, 3, 4, 5, /* 0: identity */ 866 1, 2, 0, 4, 5, 3, /* 1: a */ 867 2, 0, 1, 5, 3, 4, /* 2: a^2 */ 868 3, 4, 5, 0, 1, 2, /* 3: b */ 869 4, 5, 3, 1, 2, 0, /* 4: b a */ 870 5, 3, 4, 2, 0, 1, /* 5: b a^2 */ 871 }; 872 /* a: rotate 90 about z 873 b: swap top and bottom segments 874 r: reflect */ 875 static const PetscInt tquadVerts[16 * 8] = { 876 6, 5, 4, 7, 2, 1, 0, 3, /* -8: r b a^3 */ 877 5, 4, 7, 6, 1, 0, 3, 2, /* -7: r b a^2 */ 878 4, 7, 6, 5, 0, 3, 2, 1, /* -6: r b a */ 879 7, 6, 5, 4, 3, 2, 1, 0, /* -5: r b */ 880 2, 1, 0, 3, 6, 5, 4, 7, /* -4: r a^3 */ 881 1, 0, 3, 2, 5, 4, 7, 6, /* -3: r a^2 */ 882 0, 3, 2, 1, 4, 7, 6, 5, /* -2: r a */ 883 3, 2, 1, 0, 7, 6, 5, 4, /* -1: r */ 884 0, 1, 2, 3, 4, 5, 6, 7, /* 0: identity */ 885 1, 2, 3, 0, 5, 6, 7, 4, /* 1: a */ 886 2, 3, 0, 1, 6, 7, 4, 5, /* 2: a^2 */ 887 3, 0, 1, 2, 7, 4, 5, 6, /* 3: a^3 */ 888 4, 5, 6, 7, 0, 1, 2, 3, /* 4: b */ 889 5, 6, 7, 4, 1, 2, 3, 0, /* 5: b a */ 890 6, 7, 4, 5, 2, 3, 0, 1, /* 6: b a^2 */ 891 7, 4, 5, 6, 3, 0, 1, 2, /* 7: b a^3 */ 892 }; 893 static const PetscInt pyrVerts[8 * 5] = { 894 2, 1, 0, 3, 4, /* -4: Reflect bottom face */ 895 1, 0, 3, 2, 4, /* -3: Reflect bottom face */ 896 0, 3, 2, 1, 4, /* -2: Reflect bottom face */ 897 3, 2, 1, 0, 4, /* -1: Reflect bottom face */ 898 0, 1, 2, 3, 4, /* 0: identity */ 899 1, 2, 3, 0, 4, /* 1: 90 rotation about z */ 900 2, 3, 0, 1, 4, /* 2: 180 rotation about z */ 901 3, 0, 1, 2, 4, /* 3: 270 rotation about z */ 902 }; 903 switch (ct) { 904 case DM_POLYTOPE_POINT: 905 return pntVerts; 906 case DM_POLYTOPE_SEGMENT: 907 return &segVerts[(o + 1) * 2]; 908 case DM_POLYTOPE_POINT_PRISM_TENSOR: 909 return &segVerts[(o + 1) * 2]; 910 case DM_POLYTOPE_TRIANGLE: 911 return &triVerts[(o + 3) * 3]; 912 case DM_POLYTOPE_QUADRILATERAL: 913 return &quadVerts[(o + 4) * 4]; 914 case DM_POLYTOPE_SEG_PRISM_TENSOR: 915 return &tsegVerts[(o + 2) * 4]; 916 case DM_POLYTOPE_TETRAHEDRON: 917 return &tetVerts[(o + 12) * 4]; 918 case DM_POLYTOPE_HEXAHEDRON: 919 return &hexVerts[(o + 24) * 8]; 920 case DM_POLYTOPE_TRI_PRISM: 921 return &tripVerts[(o + 6) * 6]; 922 case DM_POLYTOPE_TRI_PRISM_TENSOR: 923 return &ttriVerts[(o + 6) * 6]; 924 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 925 return &tquadVerts[(o + 8) * 8]; 926 case DM_POLYTOPE_PYRAMID: 927 return &pyrVerts[(o + 4) * 5]; 928 default: 929 return PETSC_NULLPTR; 930 } 931 } 932 933 /* This is orientation o1 acting on orientation o2 */ 934 static inline PetscInt DMPolytopeTypeComposeOrientation(DMPolytopeType ct, PetscInt o1, PetscInt o2) 935 { 936 static const PetscInt segMult[2 * 2] = {0, -1, -1, 0}; 937 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}; 938 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, 939 -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}; 940 static const PetscInt tsegMult[4 * 4] = {0, 1, -2, -1, 1, 0, -1, -2, -2, -1, 0, 1, -1, -2, 1, 0}; 941 static const PetscInt tetMult[24 * 24] = { 942 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, 943 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, 944 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, 945 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, 946 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, 947 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, 948 -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, 949 -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, 950 -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, 951 -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, 952 -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, 953 -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, 954 }; 955 static const PetscInt hexMult[48 * 48] = { 956 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, 957 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, 958 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, 959 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, 960 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, 961 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, 962 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, 963 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, 964 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, 965 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, 966 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, 967 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, 968 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, 969 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, 970 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, 971 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, 972 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, 973 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, 974 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, 975 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, 976 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, 977 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, 978 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, 979 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, 980 -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, 981 -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, 982 -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, 983 -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, 984 -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, 985 -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, 986 -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, 987 -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, 988 -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, 989 -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, 990 -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, 991 -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, 992 -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, 993 -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, 994 -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, 995 -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, 996 -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, 997 -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, 998 -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, 999 -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, 1000 -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, 1001 -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, 1002 -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, 1003 -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, 1004 }; 1005 static const PetscInt tripMult[12 * 12] = { 1006 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, 1007 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, 1008 -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, 1009 }; 1010 static const PetscInt ttriMult[12 * 12] = { 1011 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, 1012 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, 1013 -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, 1014 }; 1015 static const PetscInt tquadMult[16 * 16] = { 1016 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, 1017 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, 1018 -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, 1019 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, 1020 -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, 1021 }; 1022 static const PetscInt pyrMult[8 * 8] = { 1023 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, 1024 }; 1025 switch (ct) { 1026 case DM_POLYTOPE_POINT: 1027 return 0; 1028 case DM_POLYTOPE_SEGMENT: 1029 case DM_POLYTOPE_POINT_PRISM_TENSOR: 1030 return segMult[(o1 + 1) * 2 + o2 + 1]; 1031 case DM_POLYTOPE_TRIANGLE: 1032 return triMult[(o1 + 3) * 6 + o2 + 3]; 1033 case DM_POLYTOPE_QUADRILATERAL: 1034 return quadMult[(o1 + 4) * 8 + o2 + 4]; 1035 case DM_POLYTOPE_SEG_PRISM_TENSOR: 1036 return tsegMult[(o1 + 2) * 4 + o2 + 2]; 1037 case DM_POLYTOPE_TETRAHEDRON: 1038 return tetMult[(o1 + 12) * 24 + o2 + 12]; 1039 case DM_POLYTOPE_HEXAHEDRON: 1040 return hexMult[(o1 + 24) * 48 + o2 + 24]; 1041 case DM_POLYTOPE_TRI_PRISM: 1042 return tripMult[(o1 + 6) * 12 + o2 + 6]; 1043 case DM_POLYTOPE_TRI_PRISM_TENSOR: 1044 return ttriMult[(o1 + 6) * 12 + o2 + 6]; 1045 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 1046 return tquadMult[(o1 + 8) * 16 + o2 + 8]; 1047 case DM_POLYTOPE_PYRAMID: 1048 return pyrMult[(o1 + 4) * 8 + o2 + 4]; 1049 default: 1050 return 0; 1051 } 1052 } 1053 1054 /* This is orientation o1 acting on orientation o2^{-1} */ 1055 static inline PetscInt DMPolytopeTypeComposeOrientationInv(DMPolytopeType ct, PetscInt o1, PetscInt o2) 1056 { 1057 static const PetscInt triInv[6] = {-3, -2, -1, 0, 2, 1}; 1058 static const PetscInt quadInv[8] = {-4, -3, -2, -1, 0, 3, 2, 1}; 1059 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}; 1060 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}; 1061 static const PetscInt tripInv[12] = {-5, -6, -4, -3, -2, -1, 0, 2, 1, 3, 4, 5}; 1062 static const PetscInt ttriInv[12] = {-6, -5, -4, -3, -2, -1, 0, 2, 1, 3, 5, 4}; 1063 static const PetscInt tquadInv[16] = {-8, -7, -6, -5, -4, -3, -2, -1, 0, 3, 2, 1, 4, 7, 6, 5}; 1064 static const PetscInt pyrInv[8] = {-4, -3, -2, -1, 0, 3, 2, 1}; 1065 switch (ct) { 1066 case DM_POLYTOPE_POINT: 1067 return 0; 1068 case DM_POLYTOPE_SEGMENT: 1069 case DM_POLYTOPE_POINT_PRISM_TENSOR: 1070 return DMPolytopeTypeComposeOrientation(ct, o1, o2); 1071 case DM_POLYTOPE_TRIANGLE: 1072 return DMPolytopeTypeComposeOrientation(ct, o1, triInv[o2 + 3]); 1073 case DM_POLYTOPE_QUADRILATERAL: 1074 return DMPolytopeTypeComposeOrientation(ct, o1, quadInv[o2 + 4]); 1075 case DM_POLYTOPE_SEG_PRISM_TENSOR: 1076 return DMPolytopeTypeComposeOrientation(ct, o1, o2); 1077 case DM_POLYTOPE_TETRAHEDRON: 1078 return DMPolytopeTypeComposeOrientation(ct, o1, tetInv[o2 + 12]); 1079 case DM_POLYTOPE_HEXAHEDRON: 1080 return DMPolytopeTypeComposeOrientation(ct, o1, hexInv[o2 + 24]); 1081 case DM_POLYTOPE_TRI_PRISM: 1082 return DMPolytopeTypeComposeOrientation(ct, o1, tripInv[o2 + 6]); 1083 case DM_POLYTOPE_TRI_PRISM_TENSOR: 1084 return DMPolytopeTypeComposeOrientation(ct, o1, ttriInv[o2 + 6]); 1085 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 1086 return DMPolytopeTypeComposeOrientation(ct, o1, tquadInv[o2 + 8]); 1087 case DM_POLYTOPE_PYRAMID: 1088 return DMPolytopeTypeComposeOrientation(ct, o1, pyrInv[o2 + 4]); 1089 default: 1090 return 0; 1091 } 1092 } 1093 1094 PETSC_EXTERN PetscErrorCode DMPolytopeMatchOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *); 1095 PETSC_EXTERN PetscErrorCode DMPolytopeMatchVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *); 1096 PETSC_EXTERN PetscErrorCode DMPolytopeGetOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *); 1097 PETSC_EXTERN PetscErrorCode DMPolytopeGetVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *); 1098 PETSC_EXTERN PetscErrorCode DMPolytopeInCellTest(DMPolytopeType, const PetscReal[], PetscBool *); 1099 1100 PETSC_EXTERN PetscErrorCode PetscDTComputeFaceQuadPermutation(DMPolytopeType, PetscQuadrature, PetscInt *, IS *[]); 1101 1102 #endif 1103