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 struct {
51: /* Inputs */
52: PetscInt dimEx; /* The dimension of the extruded mesh */
53: PetscInt cdim; /* The coordinate dimension of the input mesh */
54: PetscInt cdimEx; /* The coordinate dimension of the extruded mesh */
55: PetscInt layers; /* The number of extruded layers */
56: PetscReal thickness; /* The total thickness of the extruded layers */
57: PetscInt Nth; /* The number of specified thicknesses */
58: PetscReal *thicknesses; /* The input layer thicknesses */
59: PetscBool useTensor; /* Flag to create tensor cells */
60: PetscBool useNormal; /* Use input normal instead of calculating it */
61: PetscReal normal[3]; /* Surface normal from input */
62: PetscSimplePointFn *normalFunc; /* A function returning the normal at a given point */
63: PetscBool symmetric; /* Extrude layers symmetrically about the surface */
64: PetscBool periodic; /* Connect the extruded layer periodically to the beginning */
65: /* Calculated quantities */
66: PetscReal *layerPos; /* The position of each layer relative to the original surface, along the local normal direction */
67: PetscInt *Nt; /* The array of the number of target types */
68: DMPolytopeType **target; /* The array of target types */
69: PetscInt **size; /* The array of the number of each target type */
70: PetscInt **cone; /* The array of cones for each target cell */
71: PetscInt **ornt; /* The array of orientation for each target cell */
72: } DMPlexTransform_Extrude;
74: typedef struct {
75: PetscInt debug; // Debugging level
76: PetscBool useTensor; // Flag to create tensor cells
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: } DMPlexTransform_Cohesive;
84: typedef struct {
85: PetscInt dummy;
86: } DMPlexRefine_Regular;
88: typedef struct {
89: PetscInt dummy;
90: } DMPlexRefine_ToBox;
92: typedef struct {
93: PetscInt dummy;
94: } DMPlexRefine_Alfeld;
96: typedef struct {
97: DMLabel splitPoints; /* List of edges to be bisected (1) and cells to be divided (2) */
98: PetscSection secEdgeLen; /* Section for edge length field */
99: PetscReal *edgeLen; /* Storage for edge length field */
100: } DMPlexRefine_SBR;
102: typedef struct {
103: PetscInt dummy;
104: } DMPlexRefine_1D;
106: typedef struct {
107: PetscInt n; /* The number of divisions to produce, so n = 1 gives 2 new cells */
108: PetscReal r; /* The factor increase for cell height */
109: PetscScalar *h; /* The computed cell heights, based on r */
110: PetscInt *Nt; /* The array of the number of target types */
111: DMPolytopeType **target; /* The array of target types */
112: PetscInt **size; /* The array of the number of each target type */
113: PetscInt **cone; /* The array of cones for each target cell */
114: PetscInt **ornt; /* The array of orientation for each target cell */
115: } DMPlexRefine_BL;
117: PetscErrorCode DMPlexTransformSetDimensions_Internal(DMPlexTransform, DM, DM);
118: PetscErrorCode DMPlexTransformMapCoordinatesBarycenter_Internal(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
119: PetscErrorCode DMPlexTransformGetSubcellOrientation_Regular(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *);
120: PetscErrorCode DMPlexTransformCellRefine_Regular(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt *, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]);