1 #ifndef 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 DMPLEXREFINE1D "refine_1d" 18 #define DMPLEXEXTRUDE "extrude" 19 #define DMPLEXTRANSFORMFILTER "transform_filter" 20 21 PETSC_EXTERN PetscFunctionList DMPlexTransformList; 22 PETSC_EXTERN PetscErrorCode DMPlexTransformCreate(MPI_Comm, DMPlexTransform *); 23 PETSC_EXTERN PetscErrorCode DMPlexTransformSetType(DMPlexTransform, DMPlexTransformType); 24 PETSC_EXTERN PetscErrorCode DMPlexTransformGetType(DMPlexTransform, DMPlexTransformType *); 25 PETSC_EXTERN PetscErrorCode DMPlexTransformRegister(const char[], PetscErrorCode (*)(DMPlexTransform)); 26 PETSC_EXTERN PetscErrorCode DMPlexTransformRegisterAll(void); 27 PETSC_EXTERN PetscErrorCode DMPlexTransformRegisterDestroy(void); 28 PETSC_EXTERN PetscErrorCode DMPlexTransformSetFromOptions(DMPlexTransform); 29 PETSC_EXTERN PetscErrorCode DMPlexTransformSetUp(DMPlexTransform); 30 PETSC_EXTERN PetscErrorCode DMPlexTransformView(DMPlexTransform, PetscViewer); 31 PETSC_EXTERN PetscErrorCode DMPlexTransformDestroy(DMPlexTransform *); 32 33 PETSC_EXTERN PetscErrorCode DMPlexGetTransformType(DM, DMPlexTransformType *); 34 PETSC_EXTERN PetscErrorCode DMPlexSetTransformType(DM, DMPlexTransformType); 35 36 PETSC_EXTERN PetscErrorCode DMPlexTransformGetDM(DMPlexTransform, DM *); 37 PETSC_EXTERN PetscErrorCode DMPlexTransformSetDM(DMPlexTransform, DM); 38 PETSC_EXTERN PetscErrorCode DMPlexTransformSetDimensions(DMPlexTransform, DM, DM); 39 PETSC_EXTERN PetscErrorCode DMPlexTransformGetActive(DMPlexTransform, DMLabel *); 40 PETSC_EXTERN PetscErrorCode DMPlexTransformSetActive(DMPlexTransform, DMLabel); 41 PETSC_EXTERN PetscErrorCode DMPlexTransformGetTargetPoint(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt *); 42 PETSC_EXTERN PetscErrorCode DMPlexTransformGetSourcePoint(DMPlexTransform, PetscInt, DMPolytopeType *, DMPolytopeType *, PetscInt *, PetscInt *); 43 PETSC_EXTERN PetscErrorCode DMPlexTransformCellTransform(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt *, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]); 44 PETSC_EXTERN PetscErrorCode DMPlexTransformCellTransformIdentity(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt *, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]); 45 PETSC_EXTERN PetscErrorCode DMPlexTransformGetSubcellOrientationIdentity(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *); 46 PETSC_EXTERN PetscErrorCode DMPlexTransformGetSubcellOrientation(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *); 47 PETSC_EXTERN PetscErrorCode DMPlexTransformMapCoordinates(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]); 48 PETSC_EXTERN PetscErrorCode DMPlexTransformCreateDiscLabels(DMPlexTransform, DM); 49 PETSC_EXTERN PetscErrorCode DMPlexTransformApply(DMPlexTransform, DM, DM *); 50 PETSC_EXTERN PetscErrorCode DMPlexTransformGetCone(DMPlexTransform, PetscInt, const PetscInt *[], const PetscInt *[]); 51 PETSC_EXTERN PetscErrorCode DMPlexTransformGetConeOriented(DMPlexTransform, PetscInt, PetscInt, const PetscInt *[], const PetscInt *[]); 52 PETSC_EXTERN PetscErrorCode DMPlexTransformRestoreCone(DMPlexTransform, PetscInt, const PetscInt *[], const PetscInt *[]); 53 PETSC_EXTERN PetscErrorCode DMPlexTransformGetCellVertices(DMPlexTransform, DMPolytopeType, PetscInt *, PetscScalar *[]); 54 PETSC_EXTERN PetscErrorCode DMPlexTransformGetSubcellVertices(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt *[]); 55 PETSC_EXTERN PetscErrorCode DMPlexTransformAdaptLabel(DM, Vec, DMLabel, DMLabel, DM *); 56 57 PETSC_EXTERN PetscErrorCode DMPlexRefineRegularGetAffineTransforms(DMPlexTransform, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[]); 58 PETSC_EXTERN PetscErrorCode DMPlexRefineRegularGetAffineFaceTransforms(DMPlexTransform, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[], PetscReal *[]); 59 60 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetLayers(DMPlexTransform, PetscInt *); 61 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetLayers(DMPlexTransform, PetscInt); 62 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetThickness(DMPlexTransform, PetscReal *); 63 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetThickness(DMPlexTransform, PetscReal); 64 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetTensor(DMPlexTransform, PetscBool *); 65 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetTensor(DMPlexTransform, PetscBool); 66 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetSymmetric(DMPlexTransform, PetscBool *); 67 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetSymmetric(DMPlexTransform, PetscBool); 68 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetNormal(DMPlexTransform, PetscReal[]); 69 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetNormal(DMPlexTransform, const PetscReal[]); 70 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetNormalFunction(DMPlexTransform, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *)); 71 PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetThicknesses(DMPlexTransform, PetscInt, const PetscReal[]); 72 73 #endif 74