Actual source code: dmplextransformimpl.h
1: #pragma once
3: #include <petsc/private/dmpleximpl.h>
4: #include <petscdmplextransform.h>
6: PETSC_EXTERN PetscLogEvent DMPLEXTRANSFORM_SetUp;
7: PETSC_EXTERN PetscLogEvent DMPLEXTRANSFORM_Apply;
8: PETSC_EXTERN PetscLogEvent DMPLEXTRANSFORM_SetConeSizes;
9: PETSC_EXTERN PetscLogEvent DMPLEXTRANSFORM_SetCones;
10: PETSC_EXTERN PetscLogEvent DMPLEXTRANSFORM_CreateSF;
11: PETSC_EXTERN PetscLogEvent DMPLEXTRANSFORM_CreateLabels;
12: PETSC_EXTERN PetscLogEvent DMPLEXTRANSFORM_SetCoordinates;
14: typedef struct _p_DMPlexTransformOps *DMPlexTransformOps;
15: struct _p_DMPlexTransformOps {
16: PetscErrorCode (*view)(DMPlexTransform, PetscViewer);
17: PetscErrorCode (*setfromoptions)(DMPlexTransform, PetscOptionItems);
18: PetscErrorCode (*setup)(DMPlexTransform);
19: PetscErrorCode (*destroy)(DMPlexTransform);
20: PetscErrorCode (*setdimensions)(DMPlexTransform, DM, DM);
21: PetscErrorCode (*celltransform)(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt *, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]);
22: PetscErrorCode (*ordersupports)(DMPlexTransform, DM, DM);
23: PetscErrorCode (*getsubcellorientation)(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *);
24: PetscErrorCode (*mapcoordinates)(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
25: PetscErrorCode (*check)(DMPlexTransform, DM);
26: };
28: struct _p_DMPlexTransform {
29: PETSCHEADER(struct _p_DMPlexTransformOps);
30: void *data;
32: DM dm; /* This is the DM for which the transform has been computed */
33: DMLabel active; /* If not NULL, indicates points that are participating in the transform */
34: DMLabel trType; /* If not NULL, this holds the transformation type for each point */
35: PetscBool setupcalled; /* true if setup has been called */
36: PetscInt *ctOrderOld; /* [i] = ct: An array with original cell types in depth order */
37: PetscInt *ctOrderInvOld; /* [ct] = i: An array with the ordinal numbers for each original cell type */
38: PetscInt *ctStart; /* [ct]: The number for the first cell of each polytope type in the original mesh */
39: PetscInt *ctOrderNew; /* [i] = ct: An array with produced cell types in depth order */
40: PetscInt *ctOrderInvNew; /* [ct] = i: An array with the ordinal numbers for each produced cell type */
41: PetscInt *ctStartNew; /* [ctNew]: The number for the first cell of each polytope type in the new mesh */
42: PetscInt *offset; /* [ct/rt][ctNew]: The offset from ctStartNew[ctNew] in the new point numbering of a point of type ctNew produced from an old point of type ct or refine type rt */
43: PetscInt depth; /* The depth of the transformed mesh */
44: PetscInt *depthStart; /* The starting point for each depth stratum */
45: PetscInt *depthEnd; /* The starting point for the next depth stratum */
46: PetscInt *trNv; /* The number of transformed vertices in the closure of a cell of each type */
47: PetscScalar **trVerts; /* The transformed vertex coordinates in the closure of a cell of each type */
48: PetscInt ****trSubVerts; /* The indices for vertices of subcell (rct, r) in a cell of each type */
49: PetscFE *coordFE; /* Finite element for each cell type, used for localized coordinate interpolation */
50: PetscFEGeom **refGeom; /* Geometry of the reference cell for each cell type */
51: PetscReal redFactor; // Used to scale maxCell in each direction
52: /* Label construction */
53: PetscBool labelMatchStrata; /* Flag to restrict labeled points to the same cell type as parents */
54: PetscInt labelReplicaInc; /* Multiplier to create new label values for replicas v = oldv + r * repInc */
55: };
57: typedef struct {
58: PetscInt dummy;
59: } DMPlexTransform_Filter;
61: typedef enum {
62: NORMAL_DEFAULT,
63: NORMAL_INPUT,
64: NORMAL_COMPUTE,
65: NORMAL_COMPUTE_BD
66: } PlexNormalAlg;
67: PETSC_EXTERN const char *const PlexNormalAlgs[];
69: typedef struct {
70: /* Inputs */
71: PetscInt dimEx; /* The dimension of the extruded mesh */
72: PetscInt cdim; /* The coordinate dimension of the input mesh */
73: PetscInt cdimEx; /* The coordinate dimension of the extruded mesh */
74: PetscInt layers; /* The number of extruded layers */
75: PetscReal thickness; /* The total thickness of the extruded layers */
76: PetscInt Nth; /* The number of specified thicknesses */
77: PetscReal *thicknesses; /* The input layer thicknesses */
78: PetscBool useTensor; /* Flag to create tensor cells */
79: PlexNormalAlg normalAlg; /* Algorithm to use for computing normal */
80: PetscReal normal[3]; /* Surface normal from input */
81: DM dmNormal; // DM for normal field
82: Vec vecNormal; // Normal at each vertex
83: PetscSimplePointFn *normalFunc; /* A function returning the normal at a given point */
84: PetscBool symmetric; /* Extrude layers symmetrically about the surface */
85: PetscBool periodic; /* Connect the extruded layer periodically to the beginning */
86: /* Calculated quantities */
87: PetscReal *layerPos; /* The position of each layer relative to the original surface, along the local normal direction */
88: PetscInt *Nt; /* The array of the number of target types */
89: DMPolytopeType **target; /* The array of target types */
90: PetscInt **size; /* The array of the number of each target type */
91: PetscInt **cone; /* The array of cones for each target cell */
92: PetscInt **ornt; /* The array of orientation for each target cell */
93: // Borrowed storage
94: const PetscInt *degree; // The root degree of all points in the original mesh
95: } DMPlexTransform_Extrude;
97: typedef struct {
98: PetscInt debug; // Debugging level
99: PetscBool useTensor; // Flag to create tensor cells
100: PetscReal width; // The width of a cohesive cell
101: PetscInt *Nt; // The array of the number of target types
102: DMPolytopeType **target; // The array of target types
103: PetscInt **size; // The array of the number of each target type
104: PetscInt **cone; // The array of cones for each target cell
105: PetscInt **ornt; // The array of orientation for each target cell
106: } DMPlexTransform_Cohesive;
108: typedef struct {
109: PetscInt dummy;
110: } DMPlexRefine_Regular;
112: typedef struct {
113: PetscInt dummy;
114: } DMPlexRefine_ToBox;
116: typedef struct {
117: PetscBool reflect; /* Flag to reflect the transformation */
118: } DMPlexRefine_ToSimplex;
120: typedef struct {
121: PetscInt dummy;
122: } DMPlexRefine_Alfeld;
124: typedef struct {
125: DMLabel splitPoints; /* List of edges to be bisected (1) and cells to be divided (2) */
126: PetscSection secEdgeLen; /* Section for edge length field */
127: PetscReal *edgeLen; /* Storage for edge length field */
128: } DMPlexRefine_SBR;
130: typedef struct {
131: PetscInt dummy;
132: } DMPlexRefine_1D;
134: typedef struct {
135: PetscInt n; /* The number of divisions to produce, so n = 1 gives 2 new cells */
136: PetscReal r; /* The factor increase for cell height */
137: PetscScalar *h; /* The computed cell heights, based on r */
138: PetscInt *Nt; /* The array of the number of target types */
139: DMPolytopeType **target; /* The array of target types */
140: PetscInt **size; /* The array of the number of each target type */
141: PetscInt **cone; /* The array of cones for each target cell */
142: PetscInt **ornt; /* The array of orientation for each target cell */
143: } DMPlexRefine_BL;
145: PetscErrorCode DMPlexTransformSetDimensions_Internal(DMPlexTransform, DM, DM);
146: PetscErrorCode DMPlexTransformMapCoordinatesBarycenter_Internal(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
147: PetscErrorCode DMPlexTransformGetSubcellOrientation_Regular(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *);
148: PetscErrorCode DMPlexTransformCellRefine_Regular(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt *, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]);