1 #if !defined(PETSCDMPLEXTRANSFORM_H) 2 #define PETSCDMPLEXTRANSFORM_H 3 4 #include <petscdmplex.h> 5 #include <petscdmplextransformtypes.h> 6 7 PETSC_EXTERN PetscClassId DMPLEXTRANSFORM_CLASSID; 8 9 typedef const char* DMPlexTransformType; 10 #define DMPLEXREFINEREGULAR "refine_regular" 11 #define DMPLEXREFINEALFELD "refine_alfeld" 12 #define DMPLEXREFINEPOWELLSABIN "refine_powell_sabin" 13 #define DMPLEXREFINEBOUNDARYLAYER "refine_boundary_layer" 14 #define DMPLEXREFINESBR "refine_sbr" 15 #define DMPLEXREFINETOBOX "refine_tobox" 16 #define DMPLEXREFINETOSIMPLEX "refine_tosimplex" 17 #define DMPLEXEXTRUDE "extrude" 18 #define DMPLEXTRANSFORMFILTER "transform_filter" 19 20 PETSC_EXTERN PetscFunctionList DMPlexTransformList; 21 PETSC_EXTERN PetscErrorCode DMPlexTransformCreate(MPI_Comm, DMPlexTransform *); 22 PETSC_EXTERN PetscErrorCode DMPlexTransformSetType(DMPlexTransform, DMPlexTransformType); 23 PETSC_EXTERN PetscErrorCode DMPlexTransformGetType(DMPlexTransform, DMPlexTransformType *); 24 PETSC_EXTERN PetscErrorCode DMPlexTransformRegister(const char[], PetscErrorCode (*)(DMPlexTransform)); 25 PETSC_EXTERN PetscErrorCode DMPlexTransformRegisterAll(void); 26 PETSC_EXTERN PetscErrorCode DMPlexTransformRegisterDestroy(void); 27 PETSC_EXTERN PetscErrorCode DMPlexTransformSetFromOptions(DMPlexTransform); 28 PETSC_EXTERN PetscErrorCode DMPlexTransformSetUp(DMPlexTransform); 29 PETSC_EXTERN PetscErrorCode DMPlexTransformView(DMPlexTransform, PetscViewer); 30 PETSC_EXTERN PetscErrorCode DMPlexTransformDestroy(DMPlexTransform *); 31 32 PETSC_EXTERN PetscErrorCode DMPlexGetTransformType(DM, DMPlexTransformType *); 33 PETSC_EXTERN PetscErrorCode DMPlexSetTransformType(DM, DMPlexTransformType); 34 35 PETSC_EXTERN PetscErrorCode DMPlexTransformGetDM(DMPlexTransform, DM *); 36 PETSC_EXTERN PetscErrorCode DMPlexTransformSetDM(DMPlexTransform, DM); 37 PETSC_EXTERN PetscErrorCode DMPlexTransformSetDimensions(DMPlexTransform, DM, DM); 38 PETSC_EXTERN PetscErrorCode DMPlexTransformGetActive(DMPlexTransform, DMLabel *); 39 PETSC_EXTERN PetscErrorCode DMPlexTransformSetActive(DMPlexTransform, DMLabel); 40 PETSC_EXTERN PetscErrorCode DMPlexTransformGetTargetPoint(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt *); 41 PETSC_EXTERN PetscErrorCode DMPlexTransformGetSourcePoint(DMPlexTransform, PetscInt, DMPolytopeType *, DMPolytopeType *, PetscInt *, PetscInt *); 42 PETSC_EXTERN PetscErrorCode DMPlexTransformCellTransform(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt *, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]); 43 PETSC_EXTERN PetscErrorCode DMPlexTransformCellTransformIdentity(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt *, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]); 44 PETSC_EXTERN PetscErrorCode DMPlexTransformGetSubcellOrientationIdentity(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *); 45 PETSC_EXTERN PetscErrorCode DMPlexTransformGetSubcellOrientation(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *); 46 PETSC_EXTERN PetscErrorCode DMPlexTransformMapCoordinates(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]); 47 PETSC_EXTERN PetscErrorCode DMPlexTransformCreateDiscLabels(DMPlexTransform, DM); 48 PETSC_EXTERN PetscErrorCode DMPlexTransformApply(DMPlexTransform, DM, DM *); 49 PETSC_EXTERN PetscErrorCode DMPlexTransformGetCone(DMPlexTransform, PetscInt, const PetscInt *[], const PetscInt *[]); 50 PETSC_EXTERN PetscErrorCode DMPlexTransformGetConeOriented(DMPlexTransform, PetscInt, PetscInt, const PetscInt *[], const PetscInt *[]); 51 PETSC_EXTERN PetscErrorCode DMPlexTransformRestoreCone(DMPlexTransform, PetscInt, const PetscInt *[], const PetscInt *[]); 52 PETSC_EXTERN PetscErrorCode DMPlexTransformGetCellVertices(DMPlexTransform, DMPolytopeType, PetscInt *, PetscScalar *[]); 53 PETSC_EXTERN PetscErrorCode DMPlexTransformGetSubcellVertices(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt *[]); 54 PETSC_EXTERN PetscErrorCode DMPlexTransformAdaptLabel(DM, DMLabel, DM *); 55 56 PETSC_EXTERN PetscErrorCode DMPlexRefineRegularGetAffineTransforms(DMPlexTransform, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[]); 57 PETSC_EXTERN PetscErrorCode DMPlexRefineRegularGetAffineFaceTransforms(DMPlexTransform, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[], PetscReal *[]); 58 59 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetLayers(DMPlexTransform, PetscInt *); 60 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetLayers(DMPlexTransform, PetscInt); 61 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetThickness(DMPlexTransform, PetscReal *); 62 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetThickness(DMPlexTransform, PetscReal); 63 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetTensor(DMPlexTransform, PetscBool *); 64 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetTensor(DMPlexTransform, PetscBool); 65 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetSymmetric(DMPlexTransform, PetscBool *); 66 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetSymmetric(DMPlexTransform, PetscBool); 67 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetNormal(DMPlexTransform, const PetscReal[]); 68 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetThicknesses(DMPlexTransform, PetscInt, const PetscReal[]); 69 70 #endif 71