xref: /petsc/include/petscdmplextransform.h (revision d410b0cf18e1798d3d4c14858e0c2ffdbe2fea69)
1012bc364SMatthew G. Knepley #if !defined(PETSCDMPLEXTRANSFORM_H)
2012bc364SMatthew G. Knepley #define PETSCDMPLEXTRANSFORM_H
3012bc364SMatthew G. Knepley 
4012bc364SMatthew G. Knepley #include <petscdmplex.h>
5012bc364SMatthew G. Knepley #include <petscdmplextransformtypes.h>
6012bc364SMatthew G. Knepley 
7012bc364SMatthew G. Knepley PETSC_EXTERN PetscClassId DMPLEXTRANSFORM_CLASSID;
8012bc364SMatthew G. Knepley 
9012bc364SMatthew G. Knepley typedef const char* DMPlexTransformType;
10012bc364SMatthew G. Knepley #define DMPLEXREFINEREGULAR       "refine_regular"
11012bc364SMatthew G. Knepley #define DMPLEXREFINEALFELD        "refine_alfeld"
12012bc364SMatthew G. Knepley #define DMPLEXREFINEPOWELLSABIN   "refine_powell_sabin"
13012bc364SMatthew G. Knepley #define DMPLEXREFINEBOUNDARYLAYER "refine_boundary_layer"
14012bc364SMatthew G. Knepley #define DMPLEXREFINESBR           "refine_sbr"
15012bc364SMatthew G. Knepley #define DMPLEXREFINETOBOX         "refine_tobox"
16012bc364SMatthew G. Knepley #define DMPLEXREFINETOSIMPLEX     "refine_tosimplex"
17012bc364SMatthew G. Knepley #define DMPLEXEXTRUDE             "extrude"
18012bc364SMatthew G. Knepley #define DMPLEXTRANSFORMFILTER     "transform_filter"
19012bc364SMatthew G. Knepley 
20012bc364SMatthew G. Knepley PETSC_EXTERN PetscFunctionList DMPlexTransformList;
21012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformCreate(MPI_Comm, DMPlexTransform *);
22012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformSetType(DMPlexTransform, DMPlexTransformType);
23012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformGetType(DMPlexTransform, DMPlexTransformType *);
24012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformRegister(const char[], PetscErrorCode (*)(DMPlexTransform));
25012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformRegisterAll(void);
26012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformRegisterDestroy(void);
27012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformSetFromOptions(DMPlexTransform);
28012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformSetUp(DMPlexTransform);
29012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformView(DMPlexTransform, PetscViewer);
30012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformDestroy(DMPlexTransform *);
31012bc364SMatthew G. Knepley 
32*d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexGetTransformType(DM, DMPlexTransformType *);
33*d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexSetTransformType(DM, DMPlexTransformType);
34*d410b0cfSMatthew G. Knepley 
35012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformGetDM(DMPlexTransform, DM *);
36012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformSetDM(DMPlexTransform, DM);
37*d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformSetDimensions(DMPlexTransform, DM, DM);
38012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformGetActive(DMPlexTransform, DMLabel *);
39012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformSetActive(DMPlexTransform, DMLabel);
40012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformGetTargetPoint(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt *);
41012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformGetSourcePoint(DMPlexTransform, PetscInt, DMPolytopeType *, DMPolytopeType *, PetscInt *, PetscInt *);
42012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformCellTransform(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt *, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]);
43012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformCellTransformIdentity(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt *, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]);
44012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformGetSubcellOrientationIdentity(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *);
45012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformGetSubcellOrientation(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *);
46*d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformMapCoordinates(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
47012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformCreateDiscLabels(DMPlexTransform, DM);
48012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformApply(DMPlexTransform, DM, DM *);
49012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformGetCone(DMPlexTransform, PetscInt, const PetscInt *[], const PetscInt *[]);
50012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformGetConeOriented(DMPlexTransform, PetscInt, PetscInt, const PetscInt *[], const PetscInt *[]);
51012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformRestoreCone(DMPlexTransform, PetscInt, const PetscInt *[], const PetscInt *[]);
52012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformGetCellVertices(DMPlexTransform, DMPolytopeType, PetscInt *, PetscScalar *[]);
53012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformGetSubcellVertices(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt *[]);
54012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformAdaptLabel(DM, DMLabel, DM *);
55012bc364SMatthew G. Knepley 
56012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRefineRegularGetAffineTransforms(DMPlexTransform, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[]);
57012bc364SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexRefineRegularGetAffineFaceTransforms(DMPlexTransform, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[], PetscReal *[]);
58012bc364SMatthew G. Knepley 
59*d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetLayers(DMPlexTransform, PetscInt *);
60*d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetLayers(DMPlexTransform, PetscInt);
61*d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetThickness(DMPlexTransform, PetscReal *);
62*d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetThickness(DMPlexTransform, PetscReal);
63*d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetTensor(DMPlexTransform, PetscBool *);
64*d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetTensor(DMPlexTransform, PetscBool);
65*d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetSymmetric(DMPlexTransform, PetscBool *);
66*d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetSymmetric(DMPlexTransform, PetscBool);
67*d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetNormal(DMPlexTransform, const PetscReal[]);
68*d410b0cfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetThicknesses(DMPlexTransform, PetscInt, const PetscReal[]);
69*d410b0cfSMatthew G. Knepley 
70012bc364SMatthew G. Knepley #endif
71