Actual source code: petscdm.h

  1: /*
  2:       Objects to manage the interactions between the mesh data structures and the algebraic objects
  3: */
  4: #if !defined(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>

 12: PETSC_EXTERN PetscErrorCode DMInitializePackage(void);

 14: PETSC_EXTERN PetscClassId DM_CLASSID;
 15: PETSC_EXTERN PetscClassId DMLABEL_CLASSID;

 17: #define DMLOCATEPOINT_POINT_NOT_FOUND -367

 19: /*J
 20:     DMType - String with the name of a PETSc DM

 22:    Level: beginner

 24: .seealso: DMSetType(), DM
 25: J*/
 26: typedef const char* DMType;
 27: #define DMDA        "da"
 28: #define DMCOMPOSITE "composite"
 29: #define DMSLICED    "sliced"
 30: #define DMSHELL     "shell"
 31: #define DMPLEX      "plex"
 32: #define DMREDUNDANT "redundant"
 33: #define DMPATCH     "patch"
 34: #define DMMOAB      "moab"
 35: #define DMNETWORK   "network"
 36: #define DMFOREST    "forest"
 37: #define DMP4EST     "p4est"
 38: #define DMP8EST     "p8est"
 39: #define DMSWARM     "swarm"
 40: #define DMPRODUCT   "product"
 41: #define DMSTAG      "stag"

 43: PETSC_EXTERN const char *const DMBoundaryTypes[];
 44: PETSC_EXTERN const char *const DMBoundaryConditionTypes[];
 45: PETSC_EXTERN PetscFunctionList DMList;
 46: PETSC_EXTERN PetscErrorCode DMCreate(MPI_Comm,DM*);
 47: PETSC_EXTERN PetscErrorCode DMClone(DM,DM*);
 48: PETSC_EXTERN PetscErrorCode DMSetType(DM, DMType);
 49: PETSC_EXTERN PetscErrorCode DMGetType(DM, DMType *);
 50: PETSC_EXTERN PetscErrorCode DMRegister(const char[],PetscErrorCode (*)(DM));
 51: PETSC_EXTERN PetscErrorCode DMRegisterDestroy(void);

 53: PETSC_EXTERN PetscErrorCode DMView(DM,PetscViewer);
 54: PETSC_EXTERN PetscErrorCode DMLoad(DM,PetscViewer);
 55: PETSC_EXTERN PetscErrorCode DMDestroy(DM*);
 56: PETSC_EXTERN PetscErrorCode DMCreateGlobalVector(DM,Vec*);
 57: PETSC_EXTERN PetscErrorCode DMCreateLocalVector(DM,Vec*);
 58: PETSC_EXTERN PetscErrorCode DMGetLocalVector(DM,Vec *);
 59: PETSC_EXTERN PetscErrorCode DMRestoreLocalVector(DM,Vec *);
 60: PETSC_EXTERN PetscErrorCode DMGetGlobalVector(DM,Vec *);
 61: PETSC_EXTERN PetscErrorCode DMRestoreGlobalVector(DM,Vec *);
 62: PETSC_EXTERN PetscErrorCode DMClearGlobalVectors(DM);
 63: PETSC_EXTERN PetscErrorCode DMClearLocalVectors(DM);
 64: PETSC_EXTERN PetscErrorCode DMHasNamedGlobalVector(DM,const char*,PetscBool*);
 65: PETSC_EXTERN PetscErrorCode DMGetNamedGlobalVector(DM,const char*,Vec*);
 66: PETSC_EXTERN PetscErrorCode DMRestoreNamedGlobalVector(DM,const char*,Vec*);
 67: PETSC_EXTERN PetscErrorCode DMHasNamedLocalVector(DM,const char*,PetscBool*);
 68: PETSC_EXTERN PetscErrorCode DMGetNamedLocalVector(DM,const char*,Vec*);
 69: PETSC_EXTERN PetscErrorCode DMRestoreNamedLocalVector(DM,const char*,Vec*);
 70: PETSC_EXTERN PetscErrorCode DMGetLocalToGlobalMapping(DM,ISLocalToGlobalMapping*);
 71: PETSC_EXTERN PetscErrorCode DMCreateFieldIS(DM,PetscInt*,char***,IS**);
 72: PETSC_EXTERN PetscErrorCode DMGetBlockSize(DM,PetscInt*);
 73: PETSC_EXTERN PetscErrorCode DMCreateColoring(DM,ISColoringType,ISColoring*);
 74: PETSC_EXTERN PetscErrorCode DMCreateMatrix(DM,Mat*);
 75: PETSC_EXTERN PetscErrorCode DMSetMatrixPreallocateOnly(DM,PetscBool);
 76: PETSC_EXTERN PetscErrorCode DMSetMatrixStructureOnly(DM,PetscBool);
 77: PETSC_EXTERN PetscErrorCode DMCreateInterpolation(DM,DM,Mat*,Vec*);
 78: PETSC_EXTERN PetscErrorCode DMCreateRestriction(DM,DM,Mat*);
 79: PETSC_EXTERN PetscErrorCode DMCreateInjection(DM,DM,Mat*);
 80: PETSC_EXTERN PetscErrorCode DMCreateMassMatrix(DM,DM,Mat*);
 81: PETSC_EXTERN PetscErrorCode DMGetWorkArray(DM,PetscInt,MPI_Datatype,void*);
 82: PETSC_EXTERN PetscErrorCode DMRestoreWorkArray(DM,PetscInt,MPI_Datatype,void*);
 83: PETSC_EXTERN PetscErrorCode DMRefine(DM,MPI_Comm,DM*);
 84: PETSC_EXTERN PetscErrorCode DMCoarsen(DM,MPI_Comm,DM*);
 85: PETSC_EXTERN PetscErrorCode DMGetCoarseDM(DM,DM*);
 86: PETSC_EXTERN PetscErrorCode DMSetCoarseDM(DM,DM);
 87: PETSC_EXTERN PetscErrorCode DMGetFineDM(DM,DM*);
 88: PETSC_EXTERN PetscErrorCode DMSetFineDM(DM,DM);
 89: PETSC_EXTERN PetscErrorCode DMRefineHierarchy(DM,PetscInt,DM[]);
 90: PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy(DM,PetscInt,DM[]);
 91: PETSC_EXTERN PetscErrorCode DMCoarsenHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
 92: PETSC_EXTERN PetscErrorCode DMCoarsenHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,Vec,Mat,DM,void*),void*);
 93: PETSC_EXTERN PetscErrorCode DMRefineHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
 94: PETSC_EXTERN PetscErrorCode DMRefineHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,Mat,DM,void*),void*);
 95: PETSC_EXTERN PetscErrorCode DMRestrict(DM,Mat,Vec,Mat,DM);
 96: PETSC_EXTERN PetscErrorCode DMInterpolate(DM,Mat,DM);
 97: PETSC_EXTERN PetscErrorCode DMInterpolateSolution(DM,DM,Mat,Vec,Vec);
 98: PETSC_EXTERN PetscErrorCode DMSetFromOptions(DM);
 99: PETSC_EXTERN PetscErrorCode DMViewFromOptions(DM,PetscObject,const char[]);

101: PETSC_EXTERN PetscErrorCode DMAdaptLabel(DM,DMLabel,DM*);
102: PETSC_EXTERN PetscErrorCode DMAdaptMetric(DM, Vec, DMLabel, DM *);

104: PETSC_EXTERN PetscErrorCode DMSetUp(DM);
105: PETSC_EXTERN PetscErrorCode DMCreateInterpolationScale(DM,DM,Mat,Vec*);
106: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use DMDACreateAggregates() or DMCreateRestriction() (since version 3.12)") PetscErrorCode DMCreateAggregates(DM,DM,Mat*);
107: PETSC_EXTERN PetscErrorCode DMGlobalToLocalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
108: PETSC_EXTERN PetscErrorCode DMLocalToGlobalHookAdd(DM,PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),PetscErrorCode (*)(DM,Vec,InsertMode,Vec,void*),void*);
109: PETSC_EXTERN PetscErrorCode DMGlobalToLocal(DM,Vec,InsertMode,Vec);
110: PETSC_EXTERN PetscErrorCode DMGlobalToLocalBegin(DM,Vec,InsertMode,Vec);
111: PETSC_EXTERN PetscErrorCode DMGlobalToLocalEnd(DM,Vec,InsertMode,Vec);
112: PETSC_EXTERN PetscErrorCode DMLocalToGlobal(DM,Vec,InsertMode,Vec);
113: PETSC_EXTERN PetscErrorCode DMLocalToGlobalBegin(DM,Vec,InsertMode,Vec);
114: PETSC_EXTERN PetscErrorCode DMLocalToGlobalEnd(DM,Vec,InsertMode,Vec);
115: PETSC_EXTERN PetscErrorCode DMLocalToLocalBegin(DM,Vec,InsertMode,Vec);
116: PETSC_EXTERN PetscErrorCode DMLocalToLocalEnd(DM,Vec,InsertMode,Vec);
117: PETSC_EXTERN PetscErrorCode DMConvert(DM,DMType,DM*);

119: /* Topology support */
120: PETSC_EXTERN PetscErrorCode DMGetDimension(DM,PetscInt*);
121: PETSC_EXTERN PetscErrorCode DMSetDimension(DM,PetscInt);
122: PETSC_EXTERN PetscErrorCode DMGetDimPoints(DM,PetscInt,PetscInt*,PetscInt*);
123: PETSC_EXTERN PetscErrorCode DMGetUseNatural(DM,PetscBool*);
124: PETSC_EXTERN PetscErrorCode DMSetUseNatural(DM,PetscBool);

126: /* Coordinate support */
127: PETSC_EXTERN PetscErrorCode DMGetCoordinateDM(DM,DM*);
128: PETSC_EXTERN PetscErrorCode DMSetCoordinateDM(DM,DM);
129: PETSC_EXTERN PetscErrorCode DMGetCoordinateDim(DM,PetscInt*);
130: PETSC_EXTERN PetscErrorCode DMSetCoordinateDim(DM,PetscInt);
131: PETSC_EXTERN PetscErrorCode DMGetCoordinateSection(DM,PetscSection*);
132: PETSC_EXTERN PetscErrorCode DMSetCoordinateSection(DM,PetscInt,PetscSection);
133: PETSC_EXTERN PetscErrorCode DMGetCoordinates(DM,Vec*);
134: PETSC_EXTERN PetscErrorCode DMSetCoordinates(DM,Vec);
135: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocal(DM,Vec*);
136: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalSetUp(DM);
137: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalNoncollective(DM,Vec*);
138: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalTuple(DM,IS,PetscSection*,Vec*);
139: PETSC_EXTERN PetscErrorCode DMSetCoordinatesLocal(DM,Vec);
140: PETSC_EXTERN PetscErrorCode DMLocatePoints(DM,Vec,DMPointLocationType,PetscSF*);
141: PETSC_EXTERN PetscErrorCode DMGetPeriodicity(DM,PetscBool*,const PetscReal**,const PetscReal**,const DMBoundaryType**);
142: PETSC_EXTERN PetscErrorCode DMSetPeriodicity(DM,PetscBool,const PetscReal[],const PetscReal[],const DMBoundaryType[]);
143: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinate(DM, const PetscScalar[], PetscBool, PetscScalar[]);
144: PETSC_EXTERN PetscErrorCode DMLocalizeCoordinates(DM);
145: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalized(DM,PetscBool*);
146: PETSC_EXTERN PetscErrorCode DMGetCoordinatesLocalizedLocal(DM,PetscBool*);
147: PETSC_EXTERN PetscErrorCode DMGetNeighbors(DM,PetscInt*,const PetscMPIInt**);
148: PETSC_EXTERN PetscErrorCode DMGetCoordinateField(DM,DMField*);
149: PETSC_EXTERN PetscErrorCode DMSetCoordinateField(DM,DMField);
150: PETSC_EXTERN PetscErrorCode DMGetBoundingBox(DM,PetscReal[],PetscReal[]);
151: PETSC_EXTERN PetscErrorCode DMGetLocalBoundingBox(DM,PetscReal[],PetscReal[]);
152: PETSC_EXTERN PetscErrorCode DMProjectCoordinates(DM,PetscFE);

154: /* block hook interface */
155: PETSC_EXTERN PetscErrorCode DMSubDomainHookAdd(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
156: PETSC_EXTERN PetscErrorCode DMSubDomainHookRemove(DM,PetscErrorCode (*)(DM,DM,void*),PetscErrorCode (*)(DM,VecScatter,VecScatter,DM,void*),void*);
157: PETSC_EXTERN PetscErrorCode DMSubDomainRestrict(DM,VecScatter,VecScatter,DM);

159: PETSC_EXTERN PetscErrorCode DMSetOptionsPrefix(DM,const char []);
160: PETSC_EXTERN PetscErrorCode DMAppendOptionsPrefix(DM,const char []);
161: PETSC_EXTERN PetscErrorCode DMGetOptionsPrefix(DM,const char*[]);
162: PETSC_EXTERN PetscErrorCode DMSetVecType(DM,VecType);
163: PETSC_EXTERN PetscErrorCode DMGetVecType(DM,VecType*);
164: PETSC_EXTERN PetscErrorCode DMSetMatType(DM,MatType);
165: PETSC_EXTERN PetscErrorCode DMGetMatType(DM,MatType*);
166: PETSC_EXTERN PetscErrorCode DMSetISColoringType(DM,ISColoringType);
167: PETSC_EXTERN PetscErrorCode DMGetISColoringType(DM,ISColoringType*);
168: PETSC_EXTERN PetscErrorCode DMSetApplicationContext(DM,void*);
169: PETSC_EXTERN PetscErrorCode DMSetApplicationContextDestroy(DM,PetscErrorCode (*)(void**));
170: PETSC_EXTERN PetscErrorCode DMGetApplicationContext(DM,void*);
171: PETSC_EXTERN PetscErrorCode DMSetVariableBounds(DM,PetscErrorCode (*)(DM,Vec,Vec));
172: PETSC_EXTERN PetscErrorCode DMHasVariableBounds(DM,PetscBool *);
173: PETSC_EXTERN PetscErrorCode DMHasColoring(DM,PetscBool *);
174: PETSC_EXTERN PetscErrorCode DMHasCreateRestriction(DM,PetscBool *);
175: PETSC_EXTERN PetscErrorCode DMHasCreateInjection(DM,PetscBool *);
176: PETSC_EXTERN PetscErrorCode DMComputeVariableBounds(DM,Vec,Vec);

178: PETSC_EXTERN PetscErrorCode DMCreateSubDM(DM, PetscInt, const PetscInt[], IS *, DM *);
179: PETSC_EXTERN PetscErrorCode DMCreateSuperDM(DM[], PetscInt, IS **, DM *);
180: PETSC_EXTERN PetscErrorCode DMCreateSectionSubDM(DM,PetscInt,const PetscInt[],IS*,DM*);
181: PETSC_EXTERN PetscErrorCode DMCreateSectionSuperDM(DM[],PetscInt,IS**,DM*);
182: PETSC_EXTERN PetscErrorCode DMCreateFieldDecomposition(DM,PetscInt*,char***,IS**,DM**);
183: PETSC_EXTERN PetscErrorCode DMCreateDomainDecomposition(DM,PetscInt*,char***,IS**,IS**,DM**);
184: PETSC_EXTERN PetscErrorCode DMCreateDomainDecompositionScatters(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);

186: PETSC_EXTERN PetscErrorCode DMGetRefineLevel(DM,PetscInt*);
187: PETSC_EXTERN PetscErrorCode DMSetRefineLevel(DM,PetscInt);
188: PETSC_EXTERN PetscErrorCode DMGetCoarsenLevel(DM,PetscInt*);
189: PETSC_EXTERN PetscErrorCode DMSetCoarsenLevel(DM,PetscInt);
190: PETSC_EXTERN PetscErrorCode DMFinalizePackage(void);

192: PETSC_EXTERN PetscErrorCode VecGetDM(Vec, DM*);
193: PETSC_EXTERN PetscErrorCode VecSetDM(Vec, DM);
194: PETSC_EXTERN PetscErrorCode MatGetDM(Mat, DM*);
195: PETSC_EXTERN PetscErrorCode MatSetDM(Mat, DM);
196: PETSC_EXTERN PetscErrorCode MatFDColoringUseDM(Mat,MatFDColoring);

198: typedef struct NLF_DAAD* NLF;

200: #define DM_FILE_CLASSID 1211221

202: /* FEM support */
203: PETSC_EXTERN PetscErrorCode DMPrintCellVector(PetscInt, const char [], PetscInt, const PetscScalar []);
204: PETSC_EXTERN PetscErrorCode DMPrintCellMatrix(PetscInt, const char [], PetscInt, PetscInt, const PetscScalar []);
205: PETSC_EXTERN PetscErrorCode DMPrintLocalVec(DM, const char [], PetscReal, Vec);

207: PETSC_EXTERN PetscErrorCode DMSetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *));
208: PETSC_EXTERN PetscErrorCode DMGetNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *));
209: PETSC_EXTERN PetscErrorCode DMSetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (*)(DM, PetscInt, PetscInt, MatNullSpace *));
210: PETSC_EXTERN PetscErrorCode DMGetNearNullSpaceConstructor(DM, PetscInt, PetscErrorCode (**)(DM, PetscInt, PetscInt, MatNullSpace *));

212: PETSC_EXTERN PetscErrorCode DMGetSection(DM, PetscSection *); /* Use DMGetLocalSection() in new code (since v3.12) */
213: PETSC_EXTERN PetscErrorCode DMSetSection(DM, PetscSection);   /* Use DMSetLocalSection() in new code (since v3.12) */
214: PETSC_EXTERN PetscErrorCode DMGetLocalSection(DM, PetscSection *);
215: PETSC_EXTERN PetscErrorCode DMSetLocalSection(DM, PetscSection);
216: PETSC_EXTERN PetscErrorCode DMGetGlobalSection(DM, PetscSection *);
217: PETSC_EXTERN PetscErrorCode DMSetGlobalSection(DM, PetscSection);
218: PETSC_STATIC_INLINE PETSC_DEPRECATED_FUNCTION("Use DMGetSection() (since v3.9)") PetscErrorCode DMGetDefaultSection(DM dm, PetscSection *s) {return DMGetSection(dm,s);}
219: PETSC_STATIC_INLINE PETSC_DEPRECATED_FUNCTION("Use DMSetSection() (since v3.9)") PetscErrorCode DMSetDefaultSection(DM dm, PetscSection s) {return DMSetSection(dm,s);}
220: PETSC_STATIC_INLINE PETSC_DEPRECATED_FUNCTION("Use DMGetGlobalSection() (since v3.9)") PetscErrorCode DMGetDefaultGlobalSection(DM dm, PetscSection *s) {return DMGetGlobalSection(dm,s);}
221: PETSC_STATIC_INLINE PETSC_DEPRECATED_FUNCTION("Use DMSetGlobalSection() (since v3.9)") PetscErrorCode DMSetDefaultGlobalSection(DM dm, PetscSection s) {return DMSetGlobalSection(dm,s);}

223: PETSC_EXTERN PetscErrorCode DMGetSectionSF(DM, PetscSF*);
224: PETSC_EXTERN PetscErrorCode DMSetSectionSF(DM, PetscSF);
225: PETSC_EXTERN PetscErrorCode DMCreateSectionSF(DM, PetscSection, PetscSection);
226: PETSC_STATIC_INLINE PETSC_DEPRECATED_FUNCTION("Use DMGetSectionSF() (since v3.12)") PetscErrorCode DMGetDefaultSF(DM dm, PetscSF *s) {return DMGetSectionSF(dm,s);}
227: PETSC_STATIC_INLINE PETSC_DEPRECATED_FUNCTION("Use DMSetSectionSF() (since v3.12)") PetscErrorCode DMSetDefaultSF(DM dm, PetscSF s) {return DMSetSectionSF(dm,s);}
228: PETSC_STATIC_INLINE PETSC_DEPRECATED_FUNCTION("Use DMCreateSectionSF() (since v3.12)") PetscErrorCode DMCreateDefaultSF(DM dm, PetscSection l, PetscSection g) {return DMCreateSectionSF(dm,l,g);}
229: PETSC_EXTERN PetscErrorCode DMGetPointSF(DM, PetscSF *);
230: PETSC_EXTERN PetscErrorCode DMSetPointSF(DM, PetscSF);

232: PETSC_EXTERN PetscErrorCode DMGetDefaultConstraints(DM, PetscSection *, Mat *);
233: PETSC_EXTERN PetscErrorCode DMSetDefaultConstraints(DM, PetscSection, Mat);

235: PETSC_EXTERN PetscErrorCode DMGetOutputDM(DM, DM *);
236: PETSC_EXTERN PetscErrorCode DMGetOutputSequenceNumber(DM, PetscInt *, PetscReal *);
237: PETSC_EXTERN PetscErrorCode DMSetOutputSequenceNumber(DM, PetscInt, PetscReal);
238: PETSC_EXTERN PetscErrorCode DMOutputSequenceLoad(DM, PetscViewer, const char *, PetscInt, PetscReal *);

240: PETSC_EXTERN PetscErrorCode DMGetNumFields(DM, PetscInt *);
241: PETSC_EXTERN PetscErrorCode DMSetNumFields(DM, PetscInt);
242: PETSC_EXTERN PetscErrorCode DMGetField(DM, PetscInt, DMLabel *, PetscObject *);
243: PETSC_EXTERN PetscErrorCode DMSetField(DM, PetscInt, DMLabel, PetscObject);
244: PETSC_EXTERN PetscErrorCode DMAddField(DM, DMLabel, PetscObject);
245: PETSC_EXTERN PetscErrorCode DMSetFieldAvoidTensor(DM, PetscInt, PetscBool);
246: PETSC_EXTERN PetscErrorCode DMGetFieldAvoidTensor(DM, PetscInt, PetscBool *);
247: PETSC_EXTERN PetscErrorCode DMClearFields(DM);
248: PETSC_EXTERN PetscErrorCode DMCopyFields(DM, DM);
249: PETSC_EXTERN PetscErrorCode DMGetAdjacency(DM, PetscInt, PetscBool *, PetscBool *);
250: PETSC_EXTERN PetscErrorCode DMSetAdjacency(DM, PetscInt, PetscBool, PetscBool);
251: PETSC_EXTERN PetscErrorCode DMGetBasicAdjacency(DM, PetscBool *, PetscBool *);
252: PETSC_EXTERN PetscErrorCode DMSetBasicAdjacency(DM, PetscBool, PetscBool);

254: PETSC_EXTERN PetscErrorCode DMGetNumDS(DM, PetscInt *);
255: PETSC_EXTERN PetscErrorCode DMGetDS(DM, PetscDS *);
256: PETSC_EXTERN PetscErrorCode DMGetCellDS(DM, PetscInt, PetscDS *);
257: PETSC_EXTERN PetscErrorCode DMGetRegionDS(DM, DMLabel, IS *, PetscDS *);
258: PETSC_EXTERN PetscErrorCode DMSetRegionDS(DM, DMLabel, IS, PetscDS);
259: PETSC_EXTERN PetscErrorCode DMGetRegionNumDS(DM, PetscInt, DMLabel *, IS *, PetscDS *);
260: PETSC_EXTERN PetscErrorCode DMSetRegionNumDS(DM, PetscInt, DMLabel, IS, PetscDS);
261: PETSC_EXTERN PetscErrorCode DMFindRegionNum(DM, PetscDS, PetscInt *);
262: PETSC_EXTERN PetscErrorCode DMCreateDS(DM);
263: PETSC_EXTERN PetscErrorCode DMClearDS(DM);
264: PETSC_EXTERN PetscErrorCode DMCopyDS(DM, DM);
265: PETSC_EXTERN PetscErrorCode DMCopyDisc(DM, DM);
266: PETSC_EXTERN PetscErrorCode DMComputeExactSolution(DM, PetscReal, Vec, Vec);
267: PETSC_EXTERN PetscErrorCode DMGetNumAuxiliaryVec(DM, PetscInt *);
268: PETSC_EXTERN PetscErrorCode DMGetAuxiliaryVec(DM, DMLabel, PetscInt, Vec *);
269: PETSC_EXTERN PetscErrorCode DMSetAuxiliaryVec(DM, DMLabel, PetscInt, Vec);
270: PETSC_EXTERN PetscErrorCode DMGetAuxiliaryLabels(DM, DMLabel[], PetscInt[]);
271: PETSC_EXTERN PetscErrorCode DMCopyAuxiliaryVec(DM, DM);

273: /*MC
274:   DMInterpolationInfo - Structure for holding information about interpolation on a mesh

276:   Level: intermediate

278:   Synopsis:
279:     comm   - The communicator
280:     dim    - The spatial dimension of points
281:     nInput - The number of input points
282:     points - The input point coordinates
283:     cells  - The cell containing each point
284:     n      - The number of local points
285:     coords - The point coordinates
286:     dof    - The number of components to interpolate

288: .seealso: DMInterpolationCreate(), DMInterpolationEvaluate(), DMInterpolationAddPoints()
289: M*/
290: struct _DMInterpolationInfo {
291:   MPI_Comm   comm;
292:   PetscInt   dim;    /*1 The spatial dimension of points */
293:   PetscInt   nInput; /* The number of input points */
294:   PetscReal *points; /* The input point coordinates */
295:   PetscInt  *cells;  /* The cell containing each point */
296:   PetscInt   n;      /* The number of local points */
297:   Vec        coords; /* The point coordinates */
298:   PetscInt   dof;    /* The number of components to interpolate */
299: };
300: typedef struct _DMInterpolationInfo *DMInterpolationInfo;

302: PETSC_EXTERN PetscErrorCode DMInterpolationCreate(MPI_Comm, DMInterpolationInfo *);
303: PETSC_EXTERN PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo, PetscInt);
304: PETSC_EXTERN PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo, PetscInt *);
305: PETSC_EXTERN PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo, PetscInt);
306: PETSC_EXTERN PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo, PetscInt *);
307: PETSC_EXTERN PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo, PetscInt, PetscReal[]);
308: PETSC_EXTERN PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo, DM, PetscBool, PetscBool);
309: PETSC_EXTERN PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo, Vec *);
310: PETSC_EXTERN PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo, Vec *);
311: PETSC_EXTERN PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo, Vec *);
312: PETSC_EXTERN PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo, DM, Vec, Vec);
313: PETSC_EXTERN PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *);

315: PETSC_EXTERN PetscErrorCode DMCreateLabel(DM, const char []);
316: PETSC_EXTERN PetscErrorCode DMGetLabelValue(DM, const char[], PetscInt, PetscInt *);
317: PETSC_EXTERN PetscErrorCode DMSetLabelValue(DM, const char[], PetscInt, PetscInt);
318: PETSC_EXTERN PetscErrorCode DMClearLabelValue(DM, const char[], PetscInt, PetscInt);
319: PETSC_EXTERN PetscErrorCode DMGetLabelSize(DM, const char[], PetscInt *);
320: PETSC_EXTERN PetscErrorCode DMGetLabelIdIS(DM, const char[], IS *);
321: PETSC_EXTERN PetscErrorCode DMGetStratumSize(DM, const char [], PetscInt, PetscInt *);
322: PETSC_EXTERN PetscErrorCode DMGetStratumIS(DM, const char [], PetscInt, IS *);
323: PETSC_EXTERN PetscErrorCode DMSetStratumIS(DM, const char [], PetscInt, IS);
324: PETSC_EXTERN PetscErrorCode DMClearLabelStratum(DM, const char[], PetscInt);
325: PETSC_EXTERN PetscErrorCode DMGetLabelOutput(DM, const char[], PetscBool *);
326: PETSC_EXTERN PetscErrorCode DMSetLabelOutput(DM, const char[], PetscBool);

328: PETSC_EXTERN PetscErrorCode DMGetNumLabels(DM, PetscInt *);
329: PETSC_EXTERN PetscErrorCode DMGetLabelName(DM, PetscInt, const char **);
330: PETSC_EXTERN PetscErrorCode DMHasLabel(DM, const char [], PetscBool *);
331: PETSC_EXTERN PetscErrorCode DMGetLabel(DM, const char *, DMLabel *);
332: PETSC_EXTERN PetscErrorCode DMSetLabel(DM, DMLabel);
333: PETSC_EXTERN PetscErrorCode DMGetLabelByNum(DM, PetscInt, DMLabel *);
334: PETSC_EXTERN PetscErrorCode DMAddLabel(DM, DMLabel);
335: PETSC_EXTERN PetscErrorCode DMRemoveLabel(DM, const char [], DMLabel *);
336: PETSC_EXTERN PetscErrorCode DMRemoveLabelBySelf(DM, DMLabel *, PetscBool);
337: PETSC_EXTERN PetscErrorCode DMCopyLabels(DM, DM, PetscCopyMode, PetscBool);

339: PETSC_EXTERN PetscErrorCode DMAddBoundary(DM, DMBoundaryConditionType, const char[], DMLabel, PetscInt, const PetscInt[], PetscInt, PetscInt, const PetscInt[], void (*)(void), void (*)(void), void *, PetscInt *);
340: PETSC_EXTERN PetscErrorCode DMIsBoundaryPoint(DM, PetscInt, PetscBool *);

342: PETSC_EXTERN PetscErrorCode DMProjectFunction(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
343: PETSC_EXTERN PetscErrorCode DMProjectFunctionLocal(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void**,InsertMode,Vec);
344: PETSC_EXTERN PetscErrorCode DMProjectFunctionLabel(DM, PetscReal, DMLabel, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **, InsertMode, Vec);
345: PETSC_EXTERN PetscErrorCode DMProjectFunctionLabelLocal(DM,PetscReal,DMLabel,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,InsertMode,Vec);
346: 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);
347: 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);
348: 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);
349: PETSC_EXTERN PetscErrorCode DMComputeL2Diff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
350: PETSC_EXTERN PetscErrorCode DMComputeL2GradientDiff(DM, PetscReal, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **, Vec, const PetscReal [], PetscReal *);
351: PETSC_EXTERN PetscErrorCode DMComputeL2FieldDiff(DM,PetscReal,PetscErrorCode(**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar *,void *),void **,Vec,PetscReal *);
352: PETSC_EXTERN PetscErrorCode DMComputeError(DM, Vec, PetscReal[], Vec *);
353: PETSC_EXTERN PetscErrorCode DMHasBasisTransform(DM,PetscBool*);
354: PETSC_EXTERN PetscErrorCode DMCopyTransform(DM, DM);

356: PETSC_EXTERN PetscErrorCode DMGetCompatibility(DM,DM,PetscBool*,PetscBool*);

358: PETSC_EXTERN PetscErrorCode DMMonitorSet(DM, PetscErrorCode (*)(DM, void *), void *, PetscErrorCode (*)(void**));
359: PETSC_EXTERN PetscErrorCode DMMonitorCancel(DM);
360: PETSC_EXTERN PetscErrorCode DMMonitorSetFromOptions(DM, const char[], const char[], const char[], PetscErrorCode (*)(DM, void *), PetscErrorCode (*)(DM, PetscViewerAndFormat *), PetscBool *);
361: PETSC_EXTERN PetscErrorCode DMMonitor(DM);

363: PETSC_STATIC_INLINE PetscInt DMPolytopeTypeGetDim(DMPolytopeType ct)
364: {
365:   switch (ct) {
366:     case DM_POLYTOPE_POINT:
367:       return 0;
368:     case DM_POLYTOPE_SEGMENT:
369:     case DM_POLYTOPE_POINT_PRISM_TENSOR:
370:       return 1;
371:     case DM_POLYTOPE_TRIANGLE:
372:     case DM_POLYTOPE_QUADRILATERAL:
373:     case DM_POLYTOPE_SEG_PRISM_TENSOR:
374:       return 2;
375:     case DM_POLYTOPE_TETRAHEDRON:
376:     case DM_POLYTOPE_HEXAHEDRON:
377:     case DM_POLYTOPE_TRI_PRISM:
378:     case DM_POLYTOPE_TRI_PRISM_TENSOR:
379:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:
380:     case DM_POLYTOPE_PYRAMID:
381:       return 3;
382:     default: return -1;
383:   }
384: }

386: PETSC_STATIC_INLINE PetscInt DMPolytopeTypeGetConeSize(DMPolytopeType ct)
387: {
388:   switch (ct) {
389:     case DM_POLYTOPE_POINT:              return 0;
390:     case DM_POLYTOPE_SEGMENT:            return 2;
391:     case DM_POLYTOPE_POINT_PRISM_TENSOR: return 2;
392:     case DM_POLYTOPE_TRIANGLE:           return 3;
393:     case DM_POLYTOPE_QUADRILATERAL:      return 4;
394:     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return 4;
395:     case DM_POLYTOPE_TETRAHEDRON:        return 4;
396:     case DM_POLYTOPE_HEXAHEDRON:         return 6;
397:     case DM_POLYTOPE_TRI_PRISM:          return 5;
398:     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return 5;
399:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return 6;
400:     case DM_POLYTOPE_PYRAMID:            return 5;
401:     default: return -1;
402:   }
403: }

405: PETSC_STATIC_INLINE PetscInt DMPolytopeTypeGetNumVertices(DMPolytopeType ct)
406: {
407:   switch (ct) {
408:     case DM_POLYTOPE_POINT:              return 1;
409:     case DM_POLYTOPE_SEGMENT:            return 2;
410:     case DM_POLYTOPE_POINT_PRISM_TENSOR: return 2;
411:     case DM_POLYTOPE_TRIANGLE:           return 3;
412:     case DM_POLYTOPE_QUADRILATERAL:      return 4;
413:     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return 4;
414:     case DM_POLYTOPE_TETRAHEDRON:        return 4;
415:     case DM_POLYTOPE_HEXAHEDRON:         return 8;
416:     case DM_POLYTOPE_TRI_PRISM:          return 6;
417:     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return 6;
418:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return 8;
419:     case DM_POLYTOPE_PYRAMID:            return 5;
420:     default: return -1;
421:   }
422: }

424: PETSC_STATIC_INLINE DMPolytopeType DMPolytopeTypeSimpleShape(PetscInt dim, PetscBool simplex)
425: {
426:   return dim == 0 ? DM_POLYTOPE_POINT :
427:         (dim == 1 ? DM_POLYTOPE_SEGMENT :
428:         (dim == 2 ? (simplex ? DM_POLYTOPE_TRIANGLE : DM_POLYTOPE_QUADRILATERAL) :
429:         (dim == 3 ? (simplex ? DM_POLYTOPE_TETRAHEDRON : DM_POLYTOPE_HEXAHEDRON) : DM_POLYTOPE_UNKNOWN)));
430: }

432: PETSC_STATIC_INLINE PetscInt DMPolytopeTypeGetNumArrangments(DMPolytopeType ct)
433: {
434:   switch (ct) {
435:     case DM_POLYTOPE_POINT:              return 1;
436:     case DM_POLYTOPE_SEGMENT:            return 2;
437:     case DM_POLYTOPE_POINT_PRISM_TENSOR: return 2;
438:     case DM_POLYTOPE_TRIANGLE:           return 6;
439:     case DM_POLYTOPE_QUADRILATERAL:      return 8;
440:     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return 4;
441:     case DM_POLYTOPE_TETRAHEDRON:        return 24;
442:     case DM_POLYTOPE_HEXAHEDRON:         return 48;
443:     case DM_POLYTOPE_TRI_PRISM:          return 12;
444:     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return 12;
445:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return 16;
446:     case DM_POLYTOPE_PYRAMID:            return 8;
447:     default: return -1;
448:   }
449: }

451: /* An arrangement is a face order combined with an orientation for each face */
452: PETSC_STATIC_INLINE const PetscInt *DMPolytopeTypeGetArrangment(DMPolytopeType ct, PetscInt o)
453: {
454:   static const PetscInt pntArr[1*2] = {0, 0};
455:   /* a: swap */
456:   static const PetscInt segArr[2*2*2] = {
457:     1, 0,  0, 0, /* -1: a */
458:     0, 0,  1, 0, /*  0: e */};
459:   /* a: swap first two
460:      b: swap last two */
461:   static const PetscInt triArr[6*3*2] = {
462:     0, -1,  2, -1,  1, -1, /* -3: b */
463:     2, -1,  1, -1,  0, -1, /* -2: aba */
464:     1, -1,  0, -1,  2, -1, /* -1: a */
465:     0,  0,  1,  0,  2,  0, /*  0: identity */
466:     1,  0,  2,  0,  0,  0, /*  1: ba */
467:     2,  0,  0,  0,  1,  0, /*  2: ab */};
468:   /* a: forward cyclic permutation
469:      b: swap first and last pairs */
470:   static const PetscInt quadArr[8*4*2] = {
471:     1, -1,  0, -1,  3, -1,  2, -1, /* -4: b */
472:     0, -1,  3, -1,  2, -1,  1, -1, /* -3: b a^3 = a b */
473:     3, -1,  2, -1,  1, -1,  0, -1, /* -2: b a^2 = a^2 b */
474:     2, -1,  1, -1,  0, -1,  3, -1, /* -1: b a   = a^3 b */
475:     0,  0,  1,  0,  2,  0,  3,  0, /*  0: identity */
476:     1,  0,  2,  0,  3,  0,  0,  0, /*  1: a */
477:     2,  0,  3,  0,  0,  0,  1,  0, /*  2: a^2 */
478:     3,  0,  0,  0,  1,  0,  2,  0, /*  3: a^3 */};
479:   /* r: rotate 180
480:      b: swap top and bottom segments */
481:   static const PetscInt tsegArr[4*4*2] = {
482:     1, -1,  0, -1,  3, -1,  2, -1, /* -2: r b */
483:     0, -1,  1, -1,  3,  0,  2,  0, /* -1: r */
484:     0,  0,  1,  0,  2,  0,  3,  0, /*  0: identity */
485:     1,  0,  0,  0,  2, -1,  3, -1, /*  1: b */};
486:   /* https://en.wikiversity.org/wiki/Symmetric_group_S4 */
487:   static const PetscInt tetArr[24*4*2] = {
488:     3, -2,  2, -3,  0, -1,  1, -1, /* -12: (1324)   p22 */
489:     3, -1,  1, -3,  2, -1,  0, -1, /* -11: (14)     p21 */
490:     3, -3,  0, -3,  1, -1,  2, -1, /* -10: (1234)   p18 */
491:     2, -1,  3, -1,  1, -3,  0, -2, /*  -9: (1423)   p17 */
492:     2, -3,  0, -1,  3, -2,  1, -3, /*  -8: (1342)   p13 */
493:     2, -2,  1, -2,  0, -2,  3, -2, /*  -7: (24)     p14 */
494:     1, -2,  0, -2,  2, -2,  3, -1, /*  -6: (34)     p6  */
495:     1, -1,  3, -3,  0, -3,  2, -2, /*  -5: (1243)   p10 */
496:     1, -3,  2, -1,  3, -1,  0, -3, /*  -4: (1432)   p9  */
497:     0, -3,  1, -1,  3, -3,  2, -3, /*  -3: (12)     p1  */
498:     0, -2,  2, -2,  1, -2,  3, -3, /*  -2: (23)     p2  */
499:     0, -1,  3, -2,  2, -3,  1, -2, /*  -1: (13)     p5  */
500:     0,  0,  1,  0,  2,  0,  3,  0, /*   0: ()       p0  */
501:     0,  1,  3,  1,  1,  2,  2,  0, /*   1: (123)    p4  */
502:     0,  2,  2,  1,  3,  0,  1,  2, /*   2: (132)    p3  */
503:     1,  2,  0,  1,  3,  1,  2,  2, /*   3: (12)(34) p7  */
504:     1,  0,  2,  0,  0,  0,  3,  1, /*   4: (243)    p8  */
505:     1,  1,  3,  2,  2,  2,  0,  0, /*   5: (143)    p11 */
506:     2,  1,  3,  0,  0,  2,  1,  0, /*   6: (13)(24) p16 */
507:     2,  2,  1,  1,  3,  2,  0,  2, /*   7: (142)    p15 */
508:     2,  0,  0,  0,  1,  0,  3,  2, /*   8: (234)    p12 */
509:     3,  2,  2,  2,  1,  1,  0,  1, /*   9: (14)(23) p23 */
510:     3,  0,  0,  2,  2,  1,  1,  1, /*  10: (134)    p19 */
511:     3,  1,  1,  2,  0,  1,  2,  1  /*  11: (124)    p20 */};
512:   /* Each rotation determines a permutation of the four diagonals, and this defines the isomorphism with S_4 */
513:   static const PetscInt hexArr[48*6*2] = {
514:     2, -3,  3, -2,  4, -2,  5, -3,  1, -3,  0, -1, /* -24: reflect bottom and use -3 on top */
515:     4, -2,  5, -2,  0, -1,  1, -4,  3, -2,  2, -3, /* -23: reflect bottom and use -3 on top */
516:     5, -3,  4, -1,  1, -2,  0, -3,  3, -4,  2, -1, /* -22: reflect bottom and use -3 on top */
517:     3, -1,  2, -4,  4, -4,  5, -1,  0, -4,  1, -4, /* -21: reflect bottom and use -3 on top */
518:     3, -3,  2, -2,  5, -1,  4, -4,  1, -1,  0, -3, /* -20: reflect bottom and use -3 on top */
519:     4, -4,  5, -4,  1, -4,  0, -1,  2, -4,  3, -1, /* -19: reflect bottom and use -3 on top */
520:     2, -1,  3, -4,  5, -3,  4, -2,  0, -2,  1, -2, /* -18: reflect bottom and use -3 on top */
521:     5, -1,  4, -3,  0, -3,  1, -2,  2, -2,  3, -3, /* -17: reflect bottom and use -3 on top */
522:     4, -3,  5, -1,  3, -2,  2, -4,  1, -4,  0, -4, /* -16: reflect bottom and use -3 on top */
523:     5, -4,  4, -4,  3, -4,  2, -2,  0, -3,  1, -1, /* -15: reflect bottom and use -3 on top */
524:     3, -4,  2, -1,  1, -1,  0, -4,  4, -4,  5, -4, /* -14: reflect bottom and use -3 on top */
525:     2, -2,  3, -3,  0, -2,  1, -3,  4, -2,  5, -2, /* -13: reflect bottom and use -3 on top */
526:     1, -3,  0, -1,  4, -1,  5, -4,  3, -1,  2, -4, /* -12: reflect bottom and use -3 on top */
527:     1, -1,  0, -3,  5, -4,  4, -1,  2, -1,  3, -4, /* -11: reflect bottom and use -3 on top */
528:     5, -2,  4, -2,  2, -2,  3, -4,  1, -2,  0, -2, /* -10: reflect bottom and use -3 on top */
529:     1, -2,  0, -2,  2, -1,  3, -1,  4, -1,  5, -3, /*  -9: reflect bottom and use -3 on top */
530:     4, -1,  5, -3,  2, -4,  3, -2,  0, -1,  1, -3, /*  -8: reflect bottom and use -3 on top */
531:     3, -2,  2, -3,  0, -4,  1, -1,  5, -1,  4, -3, /*  -7: reflect bottom and use -3 on top */
532:     1, -4,  0, -4,  3, -1,  2, -1,  5, -4,  4, -4, /*  -6: reflect bottom and use -3 on top */
533:     2, -4,  3, -1,  1, -3,  0, -2,  5, -3,  4, -1, /*  -5: reflect bottom and use -3 on top */
534:     0, -4,  1, -4,  4, -3,  5, -2,  2, -3,  3, -2, /*  -4: reflect bottom and use -3 on top */
535:     0, -3,  1, -1,  3, -3,  2, -3,  4, -3,  5, -1, /*  -3: reflect bottom and use -3 on top */
536:     0, -2,  1, -2,  5, -2,  4, -3,  3, -3,  2, -2, /*  -2: reflect bottom and use -3 on top */
537:     0, -1,  1, -3,  2, -3,  3, -3,  5, -2,  4, -2, /*  -1: reflect bottom and use -3 on top */
538:     0,  0,  1,  0,  2,  0,  3,  0,  4,  0,  5,  0, /*   0: identity */
539:     0,  1,  1,  3,  5,  3,  4,  0,  2,  0,  3,  1, /*   1: 90  rotation about z */
540:     0,  2,  1,  2,  3,  0,  2,  0,  5,  3,  4,  1, /*   2: 180 rotation about z */
541:     0,  3,  1,  1,  4,  0,  5,  3,  3,  0,  2,  1, /*   3: 270 rotation about z */
542:     2,  3,  3,  2,  1,  0,  0,  3,  4,  3,  5,  1, /*   4: 90  rotation about x */
543:     1,  3,  0,  1,  3,  2,  2,  2,  4,  2,  5,  2, /*   5: 180 rotation about x */
544:     3,  1,  2,  0,  0,  1,  1,  2,  4,  1,  5,  3, /*   6: 270 rotation about x */
545:     4,  0,  5,  0,  2,  1,  3,  3,  1,  1,  0,  3, /*   7: 90  rotation about y */
546:     1,  1,  0,  3,  2,  2,  3,  2,  5,  1,  4,  3, /*   8: 180 rotation about y */
547:     5,  1,  4,  3,  2,  3,  3,  1,  0,  0,  1,  0, /*   9: 270 rotation about y */
548:     1,  0,  0,  0,  5,  1,  4,  2,  3,  2,  2,  3, /*  10: 180 rotation about x+y */
549:     1,  2,  0,  2,  4,  2,  5,  1,  2,  2,  3,  3, /*  11: 180 rotation about x-y */
550:     2,  1,  3,  0,  0,  3,  1,  0,  5,  0,  4,  0, /*  12: 180 rotation about y+z */
551:     3,  3,  2,  2,  1,  2,  0,  1,  5,  2,  4,  2, /*  13: 180 rotation about y-z */
552:     5,  3,  4,  1,  3,  1,  2,  3,  1,  3,  0,  1, /*  14: 180 rotation about z+x */
553:     4,  2,  5,  2,  3,  3,  2,  1,  0,  2,  1,  2, /*  15: 180 rotation about z-x */
554:     5,  0,  4,  0,  0,  0,  1,  3,  3,  1,  2,  0, /*  16: 120 rotation about x+y+z (v0v6) */
555:     2,  0,  3,  1,  5,  0,  4,  3,  1,  0,  0,  0, /*  17: 240 rotation about x+y+z (v0v6) */
556:     4,  3,  5,  1,  1,  1,  0,  2,  3,  3,  2,  2, /*  18: 120 rotation about x+y-z (v4v2) */
557:     3,  2,  2,  3,  5,  2,  4,  1,  0,  1,  1,  3, /*  19: 240 rotation about x+y-z (v4v2) */
558:     3,  0,  2,  1,  4,  1,  5,  2,  1,  2,  0,  2, /*  20: 120 rotation about x-y+z (v1v5) */
559:     5,  2,  4,  2,  1,  3,  0,  0,  2,  3,  3,  2, /*  21: 240 rotation about x-y+z (v1v5) */
560:     4,  1,  5,  3,  0,  2,  1,  1,  2,  1,  3,  0, /*  22: 120 rotation about x-y-z (v7v3) */
561:     2,  2,  3,  3,  4,  3,  5,  0,  0,  3,  1,  1, /*  23: 240 rotation about x-y-z (v7v3) */
562:   };
563:   static const PetscInt tripArr[12*5*2] = {
564:     1, -3,  0, -1,  3, -1,  4, -1,  2, -1, /* -6: reflect bottom and top */
565:     1, -1,  0, -3,  4, -1,  2, -1,  3, -1, /* -5: reflect bottom and top */
566:     1, -2,  0, -2,  2, -1,  3, -1,  4, -1, /* -4: reflect bottom and top */
567:     0, -3,  1, -1,  3, -3,  2, -3,  4, -3, /* -3: reflect bottom and top */
568:     0, -2,  1, -2,  4, -3,  3, -3,  2, -3, /* -2: reflect bottom and top */
569:     0, -1,  1, -3,  2, -3,  4, -3,  3, -3, /* -1: reflect bottom and top */
570:     0,  0,  1,  0,  2,  0,  3,  0,  4,  0, /*  0: identity */
571:     0,  1,  1,  2,  4,  0,  2,  0,  3,  0, /*  1: 120 rotation about z */
572:     0,  2,  1,  1,  3,  0,  4,  0,  2,  0, /*  2: 240 rotation about z */
573:     1,  1,  0,  2,  2,  2,  4,  2,  3,  2, /*  3: 180 rotation about y of 0 */
574:     1,  0,  0,  0,  4,  2,  3,  2,  2,  2, /*  4: 180 rotation about y of 1 */
575:     1,  2,  0,  1,  3,  2,  2,  2,  4,  2, /*  5: 180 rotation about y of 2 */
576:   };
577:   /* a: rotate 120 about z
578:      b: swap top and bottom segments
579:      r: reflect */
580:   static const PetscInt ttriArr[12*5*2] = {
581:     1, -3,  0, -3,  2, -2,  4, -2,  3, -2, /* -6: r b a^2 */
582:     1, -2,  0, -2,  4, -2,  3, -2,  2, -2, /* -5: r b a */
583:     1, -1,  0, -1,  3, -2,  2, -2,  4, -2, /* -4: r b */
584:     0, -3,  1, -3,  2, -1,  4, -1,  3, -1, /* -3: r a^2 */
585:     0, -2,  1, -2,  4, -1,  3, -1,  2, -1, /* -2: r a */
586:     0, -1,  1, -1,  3, -1,  2, -1,  4, -1, /* -1: r */
587:     0,  0,  1,  0,  2,  0,  3,  0,  4,  0, /*  0: identity */
588:     0,  1,  1,  1,  3,  0,  4,  0,  2,  0, /*  1: a */
589:     0,  2,  1,  2,  4,  0,  2,  0,  3,  0, /*  2: a^2 */
590:     1,  0,  0,  0,  2,  1,  3,  1,  4,  1, /*  3: b */
591:     1,  1,  0,  1,  3,  1,  4,  1,  2,  1, /*  4: b a */
592:     1,  2,  0,  2,  4,  1,  2,  1,  3,  1, /*  5: b a^2 */
593:   };
594:   /* a: rotate 90 about z
595:      b: swap top and bottom segments
596:      r: reflect */
597:   static const PetscInt tquadArr[16*6*2] = {
598:     1, -4,  0, -4,  3, -2,  2, -2,  5, -2,  4, -2, /* -8: r b a^3 */
599:     1, -3,  0, -3,  2, -2,  5, -2,  4, -2,  3, -2, /* -7: r b a^2 */
600:     1, -2,  0, -2,  5, -2,  4, -2,  3, -2,  2, -2, /* -6: r b a */
601:     1, -1,  0, -1,  4, -2,  3, -2,  2, -2,  5, -2, /* -5: r b */
602:     0, -4,  1, -4,  3, -1,  2, -1,  5, -1,  4, -1, /* -4: r a^3 */
603:     0, -3,  1, -3,  2, -1,  5, -1,  4, -1,  3, -1, /* -3: r a^2 */
604:     0, -2,  1, -2,  5, -1,  4, -1,  3, -1,  2, -1, /* -2: r a */
605:     0, -1,  1, -1,  4, -1,  3, -1,  2, -1,  5, -1, /* -1: r */
606:     0,  0,  1,  0,  2,  0,  3,  0,  4,  0,  5,  0, /*  0: identity */
607:     0,  1,  1,  1,  3,  0,  4,  0,  5,  0,  2,  0, /*  1: a */
608:     0,  2,  1,  2,  4,  0,  5,  0,  2,  0,  3,  0, /*  2: a^2 */
609:     0,  3,  1,  3,  5,  0,  2,  0,  3,  0,  4,  0, /*  3: a^3 */
610:     1,  0,  0,  0,  2,  1,  3,  1,  4,  1,  5,  1, /*  4: b */
611:     1,  1,  0,  1,  3,  1,  4,  1,  5,  1,  2,  1, /*  5: b a */
612:     1,  2,  0,  2,  4,  1,  5,  1,  2,  1,  3,  1, /*  6: b a^2 */
613:     1,  3,  0,  3,  5,  1,  2,  1,  3,  1,  4,  1, /*  7: b a^3 */
614:   };
615:   static const PetscInt pyrArr[8*5*2] = {
616:     0, -4,  2, -3,  1, -3,  4, -3,  3, -3, /* -4: Reflect bottom face */
617:     0, -3,  3, -3,  2, -3,  1, -3,  4, -3, /* -3: Reflect bottom face */
618:     0, -2,  4, -3,  3, -3,  2, -3,  1, -3, /* -2: Reflect bottom face */
619:     0, -1,  1, -3,  4, -3,  3, -3,  2, -3, /* -1: Reflect bottom face */
620:     0,  0,  1,  0,  2,  0,  3,  0,  4,  0, /*  0: identity */
621:     0,  1,  4,  0,  1,  0,  2,  0,  3,  0, /*  1:  90 rotation about z */
622:     0,  2,  3,  0,  4,  0,  1,  0,  2,  0, /*  2: 180 rotation about z */
623:     0,  3,  2,  0,  3,  0,  4,  0,  1,  0, /*  3: 270 rotation about z */
624:   };
625:   switch (ct) {
626:     case DM_POLYTOPE_POINT:              return pntArr;
627:     case DM_POLYTOPE_SEGMENT:            return &segArr[(o+1)*2*2];
628:     case DM_POLYTOPE_POINT_PRISM_TENSOR: return &segArr[(o+1)*2*2];
629:     case DM_POLYTOPE_TRIANGLE:           return &triArr[(o+3)*3*2];
630:     case DM_POLYTOPE_QUADRILATERAL:      return &quadArr[(o+4)*4*2];
631:     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return &tsegArr[(o+2)*4*2];
632:     case DM_POLYTOPE_TETRAHEDRON:        return &tetArr[(o+12)*4*2];
633:     case DM_POLYTOPE_HEXAHEDRON:         return &hexArr[(o+24)*6*2];
634:     case DM_POLYTOPE_TRI_PRISM:          return &tripArr[(o+6)*5*2];
635:     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return &ttriArr[(o+6)*5*2];
636:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return &tquadArr[(o+8)*6*2];
637:     case DM_POLYTOPE_PYRAMID:            return &pyrArr[(o+4)*5*2];
638:     default: return NULL;
639:   }
640: }

642: /* A vertex arrangment is a vertex order */
643: PETSC_STATIC_INLINE const PetscInt *DMPolytopeTypeGetVertexArrangment(DMPolytopeType ct, PetscInt o)
644: {
645:   static const PetscInt pntVerts[1]    = {0};
646:   static const PetscInt segVerts[2*2]  = {
647:     1, 0,
648:     0, 1};
649:   static const PetscInt triVerts[6*3]  = {
650:     1, 0, 2,
651:     0, 2, 1,
652:     2, 1, 0,
653:     0, 1, 2,
654:     1, 2, 0,
655:     2, 0, 1};
656:   static const PetscInt quadVerts[8*4]  = {
657:     2, 1, 0, 3,
658:     1, 0, 3, 2,
659:     0, 3, 2, 1,
660:     3, 2, 1, 0,
661:     0, 1, 2, 3,
662:     1, 2, 3, 0,
663:     2, 3, 0, 1,
664:     3, 0, 1, 2};
665:   static const PetscInt tsegVerts[4*4]  = {
666:     3, 2, 1, 0,
667:     1, 0, 3, 2,
668:     0, 1, 2, 3,
669:     2, 3, 0, 1};
670:   static const PetscInt tetVerts[24*4] = {
671:     2, 3, 1, 0, /* -12: (1324)   p22 */
672:     3, 1, 2, 0, /* -11: (14)     p21 */
673:     1, 2, 3, 0, /* -10: (1234)   p18 */
674:     3, 2, 0, 1, /*  -9: (1423)   p17 */
675:     2, 0, 3, 1, /*  -8: (1342)   p13 */
676:     0, 3, 2, 1, /*  -7: (24)     p14 */
677:     0, 1, 3, 2, /*  -6: (34)     p6  */
678:     1, 3, 0, 2, /*  -5: (1243)   p10 */
679:     3, 0, 1, 2, /*  -4: (1432    p9  */
680:     1, 0, 2, 3, /*  -3: (12)     p1  */
681:     0, 2, 1, 3, /*  -2: (23)     p2  */
682:     2, 1, 0, 3, /*  -1: (13)     p5  */
683:     0, 1, 2, 3, /*   0: ()       p0  */
684:     1, 2, 0, 3, /*   1: (123)    p4  */
685:     2, 0, 1, 3, /*   2: (132)    p3  */
686:     1, 0, 3, 2, /*   3: (12)(34) p7  */
687:     0, 3, 1, 2, /*   4: (243)    p8  */
688:     3, 1, 0, 2, /*   5: (143)    p11 */
689:     2, 3, 0, 1, /*   6: (13)(24) p16 */
690:     3, 0, 2, 1, /*   7: (142)    p15 */
691:     0, 2, 3, 1, /*   8: (234)    p12 */
692:     3, 2, 1, 0, /*   9: (14)(23) p23 */
693:     2, 1, 3, 0, /*  10: (134)    p19 */
694:     1, 3, 2, 0  /*  11: (124)    p20 */};
695:   static const PetscInt hexVerts[48*8] = {
696:     3,  0,  4,  5,  2,  6,  7,  1, /* -24: reflected 23 */
697:     3,  5,  6,  2,  0,  1,  7,  4, /* -23: reflected 22 */
698:     4,  0,  1,  7,  5,  6,  2,  3, /* -22: reflected 21 */
699:     6,  7,  1,  2,  5,  3,  0,  4, /* -21: reflected 20 */
700:     1,  2,  6,  7,  0,  4,  5,  3, /* -20: reflected 19 */
701:     6,  2,  3,  5,  7,  4,  0,  1, /* -19: reflected 18 */
702:     4,  5,  3,  0,  7,  1,  2,  6, /* -18: reflected 17 */
703:     1,  7,  4,  0,  2,  3,  5,  6, /* -17: reflected 16 */
704:     2,  3,  5,  6,  1,  7,  4,  0, /* -16: reflected 15 */
705:     7,  4,  0,  1,  6,  2,  3,  5, /* -15: reflected 14 */
706:     7,  1,  2,  6,  4,  5,  3,  0, /* -14: reflected 13 */
707:     0,  4,  5,  3,  1,  2,  6,  7, /* -13: reflected 12 */
708:     5,  4,  7,  6,  3,  2,  1,  0, /* -12: reflected 11 */
709:     7,  6,  5,  4,  1,  0,  3,  2, /* -11: reflected 10 */
710:     0,  1,  7,  4,  3,  5,  6,  2, /* -10: reflected  9 */
711:     4,  7,  6,  5,  0,  3,  2,  1, /*  -9: reflected  8 */
712:     5,  6,  2,  3,  4,  0,  1,  7, /*  -8: reflected  7 */
713:     2,  6,  7,  1,  3,  0,  4,  5, /*  -7: reflected  6 */
714:     6,  5,  4,  7,  2,  1,  0,  3, /*  -6: reflected  5 */
715:     5,  3,  0,  4,  6,  7,  1,  2, /*  -5: reflected  4 */
716:     2,  1,  0,  3,  6,  5,  4,  7, /*  -4: reflected  3 */
717:     1,  0,  3,  2,  7,  6,  5,  4, /*  -3: reflected  2 */
718:     0,  3,  2,  1,  4,  7,  6,  5, /*  -2: reflected  1 */
719:     3,  2,  1,  0,  5,  4,  7,  6, /*  -1: reflected  0 */
720:     0,  1,  2,  3,  4,  5,  6,  7, /*   0: identity */
721:     1,  2,  3,  0,  7,  4,  5,  6, /*   1: 90  rotation about z */
722:     2,  3,  0,  1,  6,  7,  4,  5, /*   2: 180 rotation about z */
723:     3,  0,  1,  2,  5,  6,  7,  4, /*   3: 270 rotation about z */
724:     4,  0,  3,  5,  7,  6,  2,  1, /*   4: 90  rotation about x */
725:     7,  4,  5,  6,  1,  2,  3,  0, /*   5: 180 rotation about x */
726:     1,  7,  6,  2,  0,  3,  5,  4, /*   6: 270 rotation about x */
727:     3,  2,  6,  5,  0,  4,  7,  1, /*   7: 90  rotation about y */
728:     5,  6,  7,  4,  3,  0,  1,  2, /*   8: 180 rotation about y */
729:     4,  7,  1,  0,  5,  3,  2,  6, /*   9: 270 rotation about y */
730:     4,  5,  6,  7,  0,  1,  2,  3, /*  10: 180 rotation about x+y */
731:     6,  7,  4,  5,  2,  3,  0,  1, /*  11: 180 rotation about x-y */
732:     3,  5,  4,  0,  2,  1,  7,  6, /*  12: 180 rotation about y+z */
733:     6,  2,  1,  7,  5,  4,  0,  3, /*  13: 180 rotation about y-z */
734:     1,  0,  4,  7,  2,  6,  5,  3, /*  14: 180 rotation about z+x */
735:     6,  5,  3,  2,  7,  1,  0,  4, /*  15: 180 rotation about z-x */
736:     0,  4,  7,  1,  3,  2,  6,  5, /*  16: 120 rotation about x+y+z (v0v6) */
737:     0,  3,  5,  4,  1,  7,  6,  2, /*  17: 240 rotation about x+y+z (v0v6) */
738:     5,  3,  2,  6,  4,  7,  1,  0, /*  18: 120 rotation about x+y-z (v4v2) */
739:     7,  6,  2,  1,  4,  0,  3,  5, /*  19: 240 rotation about x+y-z (v4v2) */
740:     2,  1,  7,  6,  3,  5,  4,  0, /*  20: 120 rotation about x-y+z (v1v5) */
741:     7,  1,  0,  4,  6,  5,  3,  2, /*  21: 240 rotation about x-y+z (v1v5) */
742:     2,  6,  5,  3,  1,  0,  4,  7, /*  22: 120 rotation about x-y-z (v7v3) */
743:     5,  4,  0,  3,  6,  2,  1,  7, /*  23: 240 rotation about x-y-z (v7v3) */
744:   };
745:   static const PetscInt tripVerts[12*6] = {
746:     4,  3,  5,  2,  1,  0, /* -6: reflect bottom and top */
747:     5,  4,  3,  1,  0,  2, /* -5: reflect bottom and top */
748:     3,  5,  4,  0,  2,  1, /* -4: reflect bottom and top */
749:     1,  0,  2,  5,  4,  3, /* -3: reflect bottom and top */
750:     0,  2,  1,  3,  5,  4, /* -2: reflect bottom and top */
751:     2,  1,  0,  4,  3,  5, /* -1: reflect bottom and top */
752:     0,  1,  2,  3,  4,  5, /*  0: identity */
753:     1,  2,  0,  5,  3,  4, /*  1: 120 rotation about z */
754:     2,  0,  1,  4,  5,  3, /*  2: 240 rotation about z */
755:     4,  5,  3,  2,  0,  1, /*  3: 180 rotation about y of 0 */
756:     3,  4,  5,  0,  1,  2, /*  4: 180 rotation about y of 1 */
757:     5,  3,  4,  1,  2,  0, /*  5: 180 rotation about y of 2 */
758:   };
759:   static const PetscInt ttriVerts[12*6] = {
760:     4,  3,  5,  1,  0,  2, /* -6: r b a^2 */
761:     3,  5,  4,  0,  2,  1, /* -5: r b a */
762:     5,  4,  3,  2,  1,  0, /* -4: r b */
763:     1,  0,  2,  4,  3,  5, /* -3: r a^2 */
764:     0,  2,  1,  3,  5,  4, /* -2: r a */
765:     2,  1,  0,  5,  4,  3, /* -1: r */
766:     0,  1,  2,  3,  4,  5, /*  0: identity */
767:     1,  2,  0,  4,  5,  3, /*  1: a */
768:     2,  0,  1,  5,  3,  4, /*  2: a^2 */
769:     3,  4,  5,  0,  1,  2, /*  3: b */
770:     4,  5,  3,  1,  2,  0, /*  4: b a */
771:     5,  3,  4,  2,  0,  1, /*  5: b a^2 */
772:   };
773:   /* a: rotate 90 about z
774:      b: swap top and bottom segments
775:      r: reflect */
776:   static const PetscInt tquadVerts[16*8] = {
777:     6,  5,  4,  7,  2,  1,  0,  3, /* -8: r b a^3 */
778:     5,  4,  7,  6,  1,  0,  3,  2, /* -7: r b a^2 */
779:     4,  7,  6,  5,  0,  3,  2,  1, /* -6: r b a */
780:     7,  6,  5,  4,  3,  2,  1,  0, /* -5: r b */
781:     2,  1,  0,  3,  6,  5,  4,  7, /* -4: r a^3 */
782:     1,  0,  3,  2,  5,  4,  7,  6, /* -3: r a^2 */
783:     0,  3,  2,  1,  4,  7,  6,  5, /* -2: r a */
784:     3,  2,  1,  0,  7,  6,  5,  4, /* -1: r */
785:     0,  1,  2,  3,  4,  5,  6,  7, /*  0: identity */
786:     1,  2,  3,  0,  5,  6,  7,  4, /*  1: a */
787:     2,  3,  0,  1,  6,  7,  4,  5, /*  2: a^2 */
788:     3,  0,  1,  2,  7,  4,  5,  6, /*  3: a^3 */
789:     4,  5,  6,  7,  0,  1,  2,  3, /*  4: b */
790:     5,  6,  7,  4,  1,  2,  3,  0, /*  5: b a */
791:     6,  7,  4,  5,  2,  3,  0,  1, /*  6: b a^2 */
792:     7,  4,  5,  6,  3,  0,  1,  2, /*  7: b a^3 */
793:   };
794:   static const PetscInt pyrVerts[8*5] = {
795:     2,  1,  0,  3,  4, /* -4: Reflect bottom face */
796:     1,  0,  3,  2,  4, /* -3: Reflect bottom face */
797:     0,  3,  2,  1,  4, /* -2: Reflect bottom face */
798:     3,  2,  1,  0,  4, /* -1: Reflect bottom face */
799:     0,  1,  2,  3,  4, /*  0: identity */
800:     1,  2,  3,  0,  4, /*  1:  90 rotation about z */
801:     2,  3,  0,  1,  4, /*  2: 180 rotation about z */
802:     3,  0,  1,  2,  4, /*  3: 270 rotation about z */
803:   };
804:   switch (ct) {
805:     case DM_POLYTOPE_POINT:              return pntVerts;
806:     case DM_POLYTOPE_SEGMENT:            return &segVerts[(o+1)*2];
807:     case DM_POLYTOPE_POINT_PRISM_TENSOR: return &segVerts[(o+1)*2];
808:     case DM_POLYTOPE_TRIANGLE:           return &triVerts[(o+3)*3];
809:     case DM_POLYTOPE_QUADRILATERAL:      return &quadVerts[(o+4)*4];
810:     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return &tsegVerts[(o+2)*4];
811:     case DM_POLYTOPE_TETRAHEDRON:        return &tetVerts[(o+12)*4];
812:     case DM_POLYTOPE_HEXAHEDRON:         return &hexVerts[(o+24)*8];
813:     case DM_POLYTOPE_TRI_PRISM:          return &tripVerts[(o+6)*6];
814:     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return &ttriVerts[(o+6)*6];
815:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return &tquadVerts[(o+8)*8];
816:     case DM_POLYTOPE_PYRAMID:            return &pyrVerts[(o+4)*5];
817:     default: return NULL;
818:   }
819: }

821: /* This is orientation o1 acting on orientation o2 */
822: PETSC_STATIC_INLINE PetscInt DMPolytopeTypeComposeOrientation(DMPolytopeType ct, PetscInt o1, PetscInt o2)
823: {
824:   static const PetscInt segMult[2*2] = {
825:      0, -1,
826:     -1,  0};
827:   static const PetscInt triMult[6*6] = {
828:      0,  2,  1, -3, -1, -2,
829:      1,  0,  2, -2, -3, -1,
830:      2,  1,  0, -1, -2, -3,
831:     -3, -2, -1,  0,  1,  2,
832:     -2, -1, -3,  1,  2,  0,
833:     -1, -3, -2,  2,  0,  1};
834:   static const PetscInt quadMult[8*8] = {
835:      0,  3,  2,  1, -4, -1, -2, -3,
836:      1,  0,  3,  2, -3, -4, -1, -2,
837:      2,  1,  0,  3, -2, -3, -4, -1,
838:      3,  2,  1,  0, -1, -2, -3, -4,
839:     -4, -3, -2, -1,  0,  1,  2,  3,
840:     -3, -2, -1, -4,  1,  2,  3,  0,
841:     -2, -1, -4, -3,  2,  3,  0,  1,
842:     -1, -4, -3, -2,  3,  0,  1,  2};
843:   static const PetscInt tsegMult[4*4] = {
844:      0,  1, -2, -1,
845:      1,  0, -1, -2,
846:     -2, -1,  0,  1,
847:     -1, -2,  1,  0};
848:   static const PetscInt tetMult[24*24] = {
849:     3, 2, 7, 0, 5, 10, 9, 8, 1, 6, 11, 4, -12, -7, -5, -9, -10, -2, -6, -1, -11, -3, -4, -8,
850:     4, 0, 8, 1, 3, 11, 10, 6, 2, 7, 9, 5, -11, -9, -4, -8, -12, -1, -5, -3, -10, -2, -6, -7,
851:     5, 1, 6, 2, 4, 9, 11, 7, 0, 8, 10, 3, -10, -8, -6, -7, -11, -3, -4, -2, -12, -1, -5, -9,
852:     0, 8, 4, 3, 11, 1, 6, 2, 10, 9, 5, 7, -9, -4, -11, -12, -1, -8, -3, -10, -5, -6, -7, -2,
853:     1, 6, 5, 4, 9, 2, 7, 0, 11, 10, 3, 8, -8, -6, -10, -11, -3, -7, -2, -12, -4, -5, -9, -1,
854:     2, 7, 3, 5, 10, 0, 8, 1, 9, 11, 4, 6, -7, -5, -12, -10, -2, -9, -1, -11, -6, -4, -8, -3,
855:     6, 5, 1, 9, 2, 4, 0, 11, 7, 3, 8, 10, -6, -10, -8, -3, -7, -11, -12, -4, -2, -9, -1, -5,
856:     7, 3, 2, 10, 0, 5, 1, 9, 8, 4, 6, 11, -5, -12, -7, -2, -9, -10, -11, -6, -1, -8, -3, -4,
857:     8, 4, 0, 11, 1, 3, 2, 10, 6, 5, 7, 9, -4, -11, -9, -1, -8, -12, -10, -5, -3, -7, -2, -6,
858:     9, 11, 10, 6, 8, 7, 3, 5, 4, 0, 2, 1, -3, -1, -2, -6, -4, -5, -9, -7, -8, -12, -10, -11,
859:     10, 9, 11, 7, 6, 8, 4, 3, 5, 1, 0, 2, -2, -3, -1, -5, -6, -4, -8, -9, -7, -11, -12, -10,
860:     11, 10, 9, 8, 7, 6, 5, 4, 3, 2, 1, 0, -1, -2, -3, -4, -5, -6, -7, -8, -9, -10, -11, -12,
861:     -12, -11, -10, -9, -8, -7, -6, -5, -4, -3, -2, -1, 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11,
862:     -11, -10, -12, -8, -7, -9, -5, -4, -6, -2, -1, -3, 1, 2, 0, 4, 5, 3, 7, 8, 6, 10, 11, 9,
863:     -10, -12, -11, -7, -9, -8, -4, -6, -5, -1, -3, -2, 2, 0, 1, 5, 3, 4, 8, 6, 7, 11, 9, 10,
864:     -9, -5, -1, -12, -2, -4, -3, -11, -7, -6, -8, -10, 3, 10, 8, 0, 7, 11, 9, 4, 2, 6, 1, 5,
865:     -8, -4, -3, -11, -1, -6, -2, -10, -9, -5, -7, -12, 4, 11, 6, 1, 8, 9, 10, 5, 0, 7, 2, 3,
866:     -7, -6, -2, -10, -3, -5, -1, -12, -8, -4, -9, -11, 5, 9, 7, 2, 6, 10, 11, 3, 1, 8, 0, 4,
867:     -3, -8, -4, -6, -11, -1, -9, -2, -10, -12, -5, -7, 6, 4, 11, 9, 1, 8, 0, 10, 5, 3, 7, 2,
868:     -2, -7, -6, -5, -10, -3, -8, -1, -12, -11, -4, -9, 7, 5, 9, 10, 2, 6, 1, 11, 3, 4, 8, 0,
869:     -1, -9, -5, -4, -12, -2, -7, -3, -11, -10, -6, -8, 8, 3, 10, 11, 0, 7, 2, 9, 4, 5, 6, 1,
870:     -6, -2, -7, -3, -5, -10, -12, -8, -1, -9, -11, -4, 9, 7, 5, 6, 10, 2, 3, 1, 11, 0, 4, 8,
871:     -5, -1, -9, -2, -4, -12, -11, -7, -3, -8, -10, -6, 10, 8, 3, 7, 11, 0, 4, 2, 9, 1, 5, 6,
872:     -4, -3, -8, -1, -6, -11, -10, -9, -2, -7, -12, -5, 11, 6, 4, 8, 9, 1, 5, 0, 10, 2, 3, 7,
873:     };
874:   static const PetscInt hexMult[48*48] = {
875:     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,
876:     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,
877:     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,
878:     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,
879:     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,
880:     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,
881:     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,
882:     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,
883:     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,
884:     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,
885:     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,
886:     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,
887:     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,
888:     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,
889:     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,
890:     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,
891:     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,
892:     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,
893:     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,
894:     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,
895:     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,
896:     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,
897:     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,
898:     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,
899:     -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,
900:     -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,
901:     -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,
902:     -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,
903:     -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,
904:     -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,
905:     -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,
906:     -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,
907:     -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,
908:     -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,
909:     -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,
910:     -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,
911:     -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,
912:     -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,
913:     -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,
914:     -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,
915:     -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,
916:     -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,
917:     -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,
918:     -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,
919:     -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,
920:     -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,
921:     -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,
922:     -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,
923:     };
924:   static const PetscInt tripMult[12*12] = {
925:     1, 0, 2, 3, 5, 4, -6, -4, -5, -2, -3, -1,
926:     0, 2, 1, 4, 3, 5, -5, -6, -4, -3, -1, -2,
927:     2, 1, 0, 5, 4, 3, -4, -5, -6, -1, -2, -3,
928:     4, 3, 5, 0, 2, 1, -3, -1, -2, -5, -6, -4,
929:     3, 5, 4, 1, 0, 2, -2, -3, -1, -6, -4, -5,
930:     5, 4, 3, 2, 1, 0, -1, -2, -3, -4, -5, -6,
931:     -6, -5, -4, -3, -2, -1, 0, 1, 2, 3, 4, 5,
932:     -4, -6, -5, -2, -1, -3, 1, 2, 0, 5, 3, 4,
933:     -5, -4, -6, -1, -3, -2, 2, 0, 1, 4, 5, 3,
934:     -3, -2, -1, -6, -5, -4, 3, 4, 5, 0, 1, 2,
935:     -1, -3, -2, -5, -4, -6, 4, 5, 3, 2, 0, 1,
936:     -2, -1, -3, -4, -6, -5, 5, 3, 4, 1, 2, 0,
937:   };
938:   static const PetscInt ttriMult[12*12] = {
939:     0, 2, 1, 3, 5, 4, -6, -4, -5, -3, -1, -2,
940:     1, 0, 2, 4, 3, 5, -5, -6, -4, -2, -3, -1,
941:     2, 1, 0, 5, 4, 3, -4, -5, -6, -1, -2, -3,
942:     3, 5, 4, 0, 2, 1, -3, -1, -2, -6, -4, -5,
943:     4, 3, 5, 1, 0, 2, -2, -3, -1, -5, -6, -4,
944:     5, 4, 3, 2, 1, 0, -1, -2, -3, -4, -5, -6,
945:     -6, -5, -4, -3, -2, -1, 0, 1, 2, 3, 4, 5,
946:     -5, -4, -6, -2, -1, -3, 1, 2, 0, 4, 5, 3,
947:     -4, -6, -5, -1, -3, -2, 2, 0, 1, 5, 3, 4,
948:     -3, -2, -1, -6, -5, -4, 3, 4, 5, 0, 1, 2,
949:     -2, -1, -3, -5, -4, -6, 4, 5, 3, 1, 2, 0,
950:     -1, -3, -2, -4, -6, -5, 5, 3, 4, 2, 0, 1,
951:   };
952:   static const PetscInt tquadMult[16*16] = {
953:     0, 3, 2, 1, 4, 7, 6, 5, -8, -5, -6, -7, -4, -1, -2, -3,
954:     1, 0, 3, 2, 5, 4, 7, 6, -7, -8, -5, -6, -3, -4, -1, -2,
955:     2, 1, 0, 3, 6, 5, 4, 7, -6, -7, -8, -5, -2, -3, -4, -1,
956:     3, 2, 1, 0, 7, 6, 5, 4, -5, -6, -7, -8, -1, -2, -3, -4,
957:     4, 7, 6, 5, 0, 3, 2, 1, -4, -1, -2, -3, -8, -5, -6, -7,
958:     5, 4, 7, 6, 1, 0, 3, 2, -3, -4, -1, -2, -7, -8, -5, -6,
959:     6, 5, 4, 7, 2, 1, 0, 3, -2, -3, -4, -1, -6, -7, -8, -5,
960:     7, 6, 5, 4, 3, 2, 1, 0, -1, -2, -3, -4, -5, -6, -7, -8,
961:     -8, -7, -6, -5, -4, -3, -2, -1, 0, 1, 2, 3, 4, 5, 6, 7,
962:     -7, -6, -5, -8, -3, -2, -1, -4, 1, 2, 3, 0, 5, 6, 7, 4,
963:     -6, -5, -8, -7, -2, -1, -4, -3, 2, 3, 0, 1, 6, 7, 4, 5,
964:     -5, -8, -7, -6, -1, -4, -3, -2, 3, 0, 1, 2, 7, 4, 5, 6,
965:     -4, -3, -2, -1, -8, -7, -6, -5, 4, 5, 6, 7, 0, 1, 2, 3,
966:     -3, -2, -1, -4, -7, -6, -5, -8, 5, 6, 7, 4, 1, 2, 3, 0,
967:     -2, -1, -4, -3, -6, -5, -8, -7, 6, 7, 4, 5, 2, 3, 0, 1,
968:     -1, -4, -3, -2, -5, -8, -7, -6, 7, 4, 5, 6, 3, 0, 1, 2,
969:   };
970:   static const PetscInt pyrMult[8*8] = {
971:     0, 3, 2, 1, -4, -1, -2, -3,
972:     1, 0, 3, 2, -3, -4, -1, -2,
973:     2, 1, 0, 3, -2, -3, -4, -1,
974:     3, 2, 1, 0, -1, -2, -3, -4,
975:     -4, -3, -2, -1, 0, 1, 2, 3,
976:     -3, -2, -1, -4, 1, 2, 3, 0,
977:     -2, -1, -4, -3, 2, 3, 0, 1,
978:     -1, -4, -3, -2, 3, 0, 1, 2,
979:   };
980:   switch (ct) {
981:     case DM_POLYTOPE_POINT:              return 0;
982:     case DM_POLYTOPE_SEGMENT:
983:     case DM_POLYTOPE_POINT_PRISM_TENSOR: return segMult[(o1+1)*2+o2+1];
984:     case DM_POLYTOPE_TRIANGLE:           return triMult[(o1+3)*6+o2+3];
985:     case DM_POLYTOPE_QUADRILATERAL:      return quadMult[(o1+4)*8+o2+4];
986:     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return tsegMult[(o1+2)*4+o2+2];
987:     case DM_POLYTOPE_TETRAHEDRON:        return tetMult[(o1+12)*24+o2+12];
988:     case DM_POLYTOPE_HEXAHEDRON:         return hexMult[(o1+24)*48+o2+24];
989:     case DM_POLYTOPE_TRI_PRISM:          return tripMult[(o1+6)*12+o2+6];
990:     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return ttriMult[(o1+6)*12+o2+6];
991:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return tquadMult[(o1+8)*16+o2+8];
992:     case DM_POLYTOPE_PYRAMID:            return pyrMult[(o1+4)*8+o2+4];
993:     default: return 0;
994:   }
995: }

997: /* This is orientation o1 acting on orientation o2^{-1} */
998: PETSC_STATIC_INLINE PetscInt DMPolytopeTypeComposeOrientationInv(DMPolytopeType ct, PetscInt o1, PetscInt o2)
999: {
1000:   static const PetscInt triInv[6]    = {-3, -2, -1, 0, 2, 1};
1001:   static const PetscInt quadInv[8]   = {-4, -3, -2, -1, 0, 3, 2, 1};
1002:   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};
1003:   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,
1004:                                           0,   3,   2,   1,   6,   5,   4,   9,   8,   7,  10,  11,  12,  13,  14, 15, 17, 16, 19, 18, 21, 20, 23, 22};
1005:   static const PetscInt tripInv[12]  = {-5, -6, -4, -3, -2, -1, 0, 2, 1, 3, 4, 5};
1006:   static const PetscInt ttriInv[12]  = {-6, -5, -4, -3, -2, -1, 0, 2, 1, 3, 5, 4};
1007:   static const PetscInt tquadInv[16] = {-8, -7, -6, -5, -4, -3, -2, -1, 0, 3, 2, 1, 4, 7, 6, 5};
1008:   static const PetscInt pyrInv[8]    = {-4, -3, -2, -1, 0, 3, 2, 1};
1009:   switch (ct) {
1010:     case DM_POLYTOPE_POINT:              return 0;
1011:     case DM_POLYTOPE_SEGMENT:
1012:     case DM_POLYTOPE_POINT_PRISM_TENSOR: return DMPolytopeTypeComposeOrientation(ct, o1, o2);
1013:     case DM_POLYTOPE_TRIANGLE:           return DMPolytopeTypeComposeOrientation(ct, o1, triInv[o2+3]);
1014:     case DM_POLYTOPE_QUADRILATERAL:      return DMPolytopeTypeComposeOrientation(ct, o1, quadInv[o2+4]);
1015:     case DM_POLYTOPE_SEG_PRISM_TENSOR:   return DMPolytopeTypeComposeOrientation(ct, o1, o2);
1016:     case DM_POLYTOPE_TETRAHEDRON:        return DMPolytopeTypeComposeOrientation(ct, o1, tetInv[o2+12]);
1017:     case DM_POLYTOPE_HEXAHEDRON:         return DMPolytopeTypeComposeOrientation(ct, o1, hexInv[o2+24]);
1018:     case DM_POLYTOPE_TRI_PRISM:          return DMPolytopeTypeComposeOrientation(ct, o1, tripInv[o2+6]);
1019:     case DM_POLYTOPE_TRI_PRISM_TENSOR:   return DMPolytopeTypeComposeOrientation(ct, o1, ttriInv[o2+6]);
1020:     case DM_POLYTOPE_QUAD_PRISM_TENSOR:  return DMPolytopeTypeComposeOrientation(ct, o1, tquadInv[o2+8]);
1021:     case DM_POLYTOPE_PYRAMID:            return DMPolytopeTypeComposeOrientation(ct, o1, pyrInv[o2+4]);
1022:     default: return 0;
1023:   }
1024: }

1026: PETSC_EXTERN PetscErrorCode DMPolytopeMatchOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *);
1027: PETSC_EXTERN PetscErrorCode DMPolytopeMatchVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *, PetscBool *);
1028: PETSC_EXTERN PetscErrorCode DMPolytopeGetOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *);
1029: PETSC_EXTERN PetscErrorCode DMPolytopeGetVertexOrientation(DMPolytopeType, const PetscInt[], const PetscInt[], PetscInt *);
1030: PETSC_EXTERN PetscErrorCode DMPolytopeInCellTest(DMPolytopeType, const PetscReal[], PetscBool *);

1032: #endif