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 (*ordersupports)(DMPlexTransform, DM, DM);
 15:   PetscErrorCode (*getsubcellorientation)(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *);
 16:   PetscErrorCode (*mapcoordinates)(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
 17: };

 19: struct _p_DMPlexTransform {
 20:   PETSCHEADER(struct _p_DMPlexTransformOps);
 21:   void *data;

 23:   DM            dm;            /* This is the DM for which the transform has been computed */
 24:   DMLabel       active;        /* If not NULL, indicates points that are participating in the transform */
 25:   DMLabel       trType;        /* If not NULL, this holds the transformation type for each point */
 26:   PetscBool     setupcalled;   /* true if setup has been called */
 27:   PetscInt     *ctOrderOld;    /* [i] = ct: An array with original cell types in depth order */
 28:   PetscInt     *ctOrderInvOld; /* [ct] = i: An array with the ordinal numbers for each original cell type */
 29:   PetscInt     *ctStart;       /* [ct]: The number for the first cell of each polytope type in the original mesh */
 30:   PetscInt     *ctOrderNew;    /* [i] = ct: An array with produced cell types in depth order */
 31:   PetscInt     *ctOrderInvNew; /* [ct] = i: An array with the ordinal numbers for each produced cell type */
 32:   PetscInt     *ctStartNew;    /* [ctNew]: The number for the first cell of each polytope type in the new mesh */
 33:   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 */
 34:   PetscInt      depth;         /* The depth of the transformed mesh */
 35:   PetscInt     *depthStart;    /* The starting point for each depth stratum */
 36:   PetscInt     *depthEnd;      /* The starting point for the next depth stratum */
 37:   PetscInt     *trNv;          /* The number of transformed vertices in the closure of a cell of each type */
 38:   PetscScalar **trVerts;       /* The transformed vertex coordinates in the closure of a cell of each type */
 39:   PetscInt  ****trSubVerts;    /* The indices for vertices of subcell (rct, r) in a cell of each type */
 40:   PetscFE      *coordFE;       /* Finite element for each cell type, used for localized coordinate interpolation */
 41:   PetscFEGeom **refGeom;       /* Geometry of the reference cell for each cell type */
 42:   /* Label construction */
 43:   PetscBool labelMatchStrata; /* Flag to restrict labeled points to the same cell type as parents */
 44:   PetscInt  labelReplicaInc;  /* Multiplier to create new label values for replicas v = oldv + r * repInc */
 45: };

 47: typedef struct {
 48:   PetscInt dummy;
 49: } DMPlexTransform_Filter;

 51: typedef enum {
 52:   NORMAL_DEFAULT,
 53:   NORMAL_INPUT,
 54:   NORMAL_COMPUTE,
 55:   NORMAL_COMPUTE_BD
 56: } PlexNormalAlg;
 57: PETSC_EXTERN const char *const PlexNormalAlgs[];

 59: typedef struct {
 60:   /* Inputs */
 61:   PetscInt            dimEx;       /* The dimension of the extruded mesh */
 62:   PetscInt            cdim;        /* The coordinate dimension of the input mesh */
 63:   PetscInt            cdimEx;      /* The coordinate dimension of the extruded mesh */
 64:   PetscInt            layers;      /* The number of extruded layers */
 65:   PetscReal           thickness;   /* The total thickness of the extruded layers */
 66:   PetscInt            Nth;         /* The number of specified thicknesses */
 67:   PetscReal          *thicknesses; /* The input layer thicknesses */
 68:   PetscBool           useTensor;   /* Flag to create tensor cells */
 69:   PlexNormalAlg       normalAlg;   /* Algorithm to use for computing normal */
 70:   PetscReal           normal[3];   /* Surface normal from input */
 71:   DM                  dmNormal;    // DM for normal field
 72:   Vec                 vecNormal;   // Normal at each vertex
 73:   PetscSimplePointFn *normalFunc;  /* A function returning the normal at a given point */
 74:   PetscBool           symmetric;   /* Extrude layers symmetrically about the surface */
 75:   PetscBool           periodic;    /* Connect the extruded layer periodically to the beginning */
 76:   /* Calculated quantities */
 77:   PetscReal       *layerPos; /* The position of each layer relative to the original surface, along the local normal direction */
 78:   PetscInt        *Nt;       /* The array of the number of target types */
 79:   DMPolytopeType **target;   /* The array of target types */
 80:   PetscInt       **size;     /* The array of the number of each target type */
 81:   PetscInt       **cone;     /* The array of cones for each target cell */
 82:   PetscInt       **ornt;     /* The array of orientation for each target cell */
 83:   // Borrowed storage
 84:   const PetscInt *degree; // The root degree of all points in the original mesh
 85: } DMPlexTransform_Extrude;

 87: typedef struct {
 88:   PetscInt         debug;     // Debugging level
 89:   PetscBool        useTensor; // Flag to create tensor cells
 90:   PetscReal        width;     // The width of a cohesive cell
 91:   PetscInt        *Nt;        // The array of the number of target types
 92:   DMPolytopeType **target;    // The array of target types
 93:   PetscInt       **size;      // The array of the number of each target type
 94:   PetscInt       **cone;      // The array of cones for each target cell
 95:   PetscInt       **ornt;      // The array of orientation for each target cell
 96: } DMPlexTransform_Cohesive;

 98: typedef struct {
 99:   PetscInt dummy;
100: } DMPlexRefine_Regular;

102: typedef struct {
103:   PetscInt dummy;
104: } DMPlexRefine_ToBox;

106: typedef struct {
107:   PetscInt dummy;
108: } DMPlexRefine_Alfeld;

110: typedef struct {
111:   DMLabel      splitPoints; /* List of edges to be bisected (1) and cells to be divided (2) */
112:   PetscSection secEdgeLen;  /* Section for edge length field */
113:   PetscReal   *edgeLen;     /* Storage for edge length field */
114: } DMPlexRefine_SBR;

116: typedef struct {
117:   PetscInt dummy;
118: } DMPlexRefine_1D;

120: typedef struct {
121:   PetscInt         n;      /* The number of divisions to produce, so n = 1 gives 2 new cells */
122:   PetscReal        r;      /* The factor increase for cell height */
123:   PetscScalar     *h;      /* The computed cell heights, based on r */
124:   PetscInt        *Nt;     /* The array of the number of target types */
125:   DMPolytopeType **target; /* The array of target types */
126:   PetscInt       **size;   /* The array of the number of each target type */
127:   PetscInt       **cone;   /* The array of cones for each target cell */
128:   PetscInt       **ornt;   /* The array of orientation for each target cell */
129: } DMPlexRefine_BL;

131: PetscErrorCode DMPlexTransformSetDimensions_Internal(DMPlexTransform, DM, DM);
132: PetscErrorCode DMPlexTransformMapCoordinatesBarycenter_Internal(DMPlexTransform, DMPolytopeType, DMPolytopeType, PetscInt, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
133: PetscErrorCode DMPlexTransformGetSubcellOrientation_Regular(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscInt *, PetscInt *);
134: PetscErrorCode DMPlexTransformCellRefine_Regular(DMPlexTransform, DMPolytopeType, PetscInt, PetscInt *, PetscInt *, DMPolytopeType *[], PetscInt *[], PetscInt *[], PetscInt *[]);