Actual source code: petscdmplextransform.h
1: #pragma once
3: #include <petscdmplex.h>
4: #include <petscdmplextransformtypes.h>
6: PETSC_EXTERN PetscClassId DMPLEXTRANSFORM_CLASSID;
8: /*J
9: DMPlexTransformType - String with the name of a PETSc DMPlexTransform type
11: Level: beginner
13: Note:
14: [](plex_transform_table) for a table of available matrix types
16: .seealso: [](plex_transform_table), [](ch_unstructured), `DMPlexTransformCreate()`, `DMPlexTransform`, `DMPlexTransformRegister()`
17: J*/
18: typedef const char *DMPlexTransformType;
19: #define DMPLEXREFINEREGULAR "refine_regular"
20: #define DMPLEXREFINEALFELD "refine_alfeld"
21: #define DMPLEXREFINEPOWELLSABIN "refine_powell_sabin"
22: #define DMPLEXREFINEBOUNDARYLAYER "refine_boundary_layer"
23: #define DMPLEXREFINESBR "refine_sbr"
24: #define DMPLEXREFINETOBOX "refine_tobox"
25: #define DMPLEXREFINETOSIMPLEX "refine_tosimplex"
26: #define DMPLEXREFINE1D "refine_1d"
27: #define DMPLEXEXTRUDE "extrude"
28: #define DMPLEXCOHESIVEEXTRUDE "cohesive_extrude"
29: #define DMPLEXTRANSFORMFILTER "transform_filter"
31: PETSC_EXTERN PetscFunctionList DMPlexTransformList;
32: PETSC_EXTERN PetscErrorCode DMPlexTransformCreate(MPI_Comm, DMPlexTransform *);
33: PETSC_EXTERN PetscErrorCode DMPlexTransformSetType(DMPlexTransform, DMPlexTransformType);
34: PETSC_EXTERN PetscErrorCode DMPlexTransformGetType(DMPlexTransform, DMPlexTransformType *);
35: PETSC_EXTERN PetscErrorCode DMPlexTransformRegister(const char[], PetscErrorCode (*)(DMPlexTransform));
36: PETSC_EXTERN PetscErrorCode DMPlexTransformRegisterAll(void);
37: PETSC_EXTERN PetscErrorCode DMPlexTransformRegisterDestroy(void);
38: PETSC_EXTERN PetscErrorCode DMPlexTransformSetFromOptions(DMPlexTransform);
39: PETSC_EXTERN PetscErrorCode DMPlexTransformSetUp(DMPlexTransform);
40: PETSC_EXTERN PetscErrorCode DMPlexTransformView(DMPlexTransform, PetscViewer);
41: PETSC_EXTERN PetscErrorCode DMPlexTransformDestroy(DMPlexTransform *);
43: PETSC_EXTERN PetscErrorCode DMPlexGetTransformType(DM, DMPlexTransformType *);
44: PETSC_EXTERN PetscErrorCode DMPlexSetTransformType(DM, DMPlexTransformType);
46: PETSC_EXTERN PetscErrorCode DMPlexTransformGetDM(DMPlexTransform, DM *);
47: PETSC_EXTERN PetscErrorCode DMPlexTransformSetDM(DMPlexTransform, DM);
48: PETSC_EXTERN PetscErrorCode DMPlexTransformSetDimensions(DMPlexTransform, DM, DM);
49: PETSC_EXTERN PetscErrorCode DMPlexTransformGetChart(DMPlexTransform, PetscInt *, PetscInt *);
50: PETSC_EXTERN PetscErrorCode DMPlexTransformGetCellType(DMPlexTransform, PetscInt, DMPolytopeType *);
51: PETSC_EXTERN PetscErrorCode DMPlexTransformGetCellTypeStratum(DMPlexTransform, DMPolytopeType, PetscInt *, PetscInt *);
52: PETSC_EXTERN PetscErrorCode DMPlexTransformGetDepth(DMPlexTransform, PetscInt *);
53: PETSC_EXTERN PetscErrorCode DMPlexTransformGetDepthStratum(DMPlexTransform, PetscInt, PetscInt *, PetscInt *);
54: PETSC_EXTERN PetscErrorCode DMPlexTransformGetActive(DMPlexTransform, DMLabel *);
55: PETSC_EXTERN PetscErrorCode DMPlexTransformSetActive(DMPlexTransform, DMLabel);
56: PETSC_EXTERN PetscErrorCode DMPlexTransformGetMatchStrata(DMPlexTransform, PetscBool *);
57: PETSC_EXTERN PetscErrorCode DMPlexTransformSetMatchStrata(DMPlexTransform, PetscBool);
59: PETSC_EXTERN PetscErrorCode DMPlexTransformGetTargetPoint(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt *);
60: PETSC_EXTERN PetscErrorCode DMPlexTransformGetSourcePoint(DMPlexTransform, PetscInt, DMPolytopeType *, DMPolytopeType *, PetscInt *, PetscInt *);
61: PETSC_EXTERN PetscErrorCode DMPlexTransformCellTransform(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt *, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]);
62: PETSC_EXTERN PetscErrorCode DMPlexTransformCellTransformIdentity(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt *, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]);
63: PETSC_EXTERN PetscErrorCode DMPlexTransformGetSubcellOrientationIdentity(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *);
64: PETSC_EXTERN PetscErrorCode DMPlexTransformGetSubcellOrientation(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *);
65: PETSC_EXTERN PetscErrorCode DMPlexTransformMapCoordinates(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
66: PETSC_EXTERN PetscErrorCode DMPlexTransformCreateDiscLabels(DMPlexTransform, DM);
67: PETSC_EXTERN PetscErrorCode DMPlexTransformApply(DMPlexTransform, DM, DM *);
68: PETSC_EXTERN PetscErrorCode DMPlexTransformGetConeSize(DMPlexTransform, PetscInt, PetscInt *);
69: PETSC_EXTERN PetscErrorCode DMPlexTransformGetCone(DMPlexTransform, PetscInt, const PetscInt *[], const PetscInt *[]);
70: PETSC_EXTERN PetscErrorCode DMPlexTransformGetConeOriented(DMPlexTransform, PetscInt, PetscInt, const PetscInt *[], const PetscInt *[]);
71: PETSC_EXTERN PetscErrorCode DMPlexTransformRestoreCone(DMPlexTransform, PetscInt, const PetscInt *[], const PetscInt *[]);
72: PETSC_EXTERN PetscErrorCode DMPlexTransformGetCellVertices(DMPlexTransform, DMPolytopeType, PetscInt *, PetscScalar *[]);
73: PETSC_EXTERN PetscErrorCode DMPlexTransformGetSubcellVertices(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt *[]);
74: PETSC_EXTERN PetscErrorCode DMPlexTransformAdaptLabel(DM, Vec, DMLabel, DMLabel, DM *);
76: PETSC_EXTERN PetscErrorCode DMPlexRefineRegularGetAffineTransforms(DMPlexTransform, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[]);
77: PETSC_EXTERN PetscErrorCode DMPlexRefineRegularGetAffineFaceTransforms(DMPlexTransform, DMPolytopeType, PetscInt *, PetscReal *[], PetscReal *[], PetscReal *[], PetscReal *[]);
79: PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetLayers(DMPlexTransform, PetscInt *);
80: PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetLayers(DMPlexTransform, PetscInt);
81: PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetThickness(DMPlexTransform, PetscReal *);
82: PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetThickness(DMPlexTransform, PetscReal);
83: PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetTensor(DMPlexTransform, PetscBool *);
84: PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetTensor(DMPlexTransform, PetscBool);
85: PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetSymmetric(DMPlexTransform, PetscBool *);
86: PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetSymmetric(DMPlexTransform, PetscBool);
87: PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetPeriodic(DMPlexTransform, PetscBool *);
88: PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetPeriodic(DMPlexTransform, PetscBool);
89: PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeGetNormal(DMPlexTransform, PetscReal[]);
90: PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetNormal(DMPlexTransform, const PetscReal[]);
91: PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetNormalFunction(DMPlexTransform, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *));
92: PETSC_EXTERN PetscErrorCode DMPlexTransformExtrudeSetThicknesses(DMPlexTransform, PetscInt, const PetscReal[]);
94: PETSC_EXTERN PetscErrorCode DMPlexTransformCohesiveExtrudeGetTensor(DMPlexTransform, PetscBool *);
95: PETSC_EXTERN PetscErrorCode DMPlexTransformCohesiveExtrudeSetTensor(DMPlexTransform, PetscBool);
96: PETSC_EXTERN PetscErrorCode DMPlexTransformCohesiveExtrudeGetWidth(DMPlexTransform, PetscReal *);
97: PETSC_EXTERN PetscErrorCode DMPlexTransformCohesiveExtrudeSetWidth(DMPlexTransform, PetscReal);
99: PETSC_EXTERN PetscErrorCode DMPlexCreateEphemeral(DMPlexTransform, const char[], DM *);