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