xref: /petsc/src/dm/impls/plex/plexextrude.c (revision 33e3ecb2fe4b52c32599116bc622b359dfa41247)
1 #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
2 #include <petscdmplextransform.h>
3 
4 /*@
5   DMPlexExtrude - Extrude a volumetric mesh from the input surface mesh
6 
7   Input Parameters:
8 + dm          - The surface mesh
9 . layers      - The number of extruded layers
10 . thickness   - The total thickness of the extruded layers, or `PETSC_DETERMINE`
11 . tensor      - Flag to create tensor produt cells
12 . symmetric   - Flag to extrude symmetrically about the surface
13 . periodic    - Flag to extrude periodically
14 . normal      - Surface normal vector, or `NULL`
15 - thicknesses - Thickness of each layer, or `NULL`
16 
17   Output Parameter:
18 . edm - The volumetric mesh
19 
20   Options Database Keys:
21 + -dm_plex_transform_extrude_thickness <t>           - The total thickness of extruded layers
22 . -dm_plex_transform_extrude_use_tensor <bool>       - Use tensor cells when extruding
23 . -dm_plex_transform_extrude_symmetric <bool>        - Extrude layers symmetrically about the surface
24 . -dm_plex_transform_extrude_periodic <bool>         - Extrude layers periodically
25 . -dm_plex_transform_extrude_normal <n0,...,nd>      - Specify the extrusion direction
26 - -dm_plex_transform_extrude_thicknesses <t0,...,tl> - Specify thickness of each layer
27 
28   Level: intermediate
29 
30   Notes:
31   Extrusion is implemented as a `DMPlexTransform`, so that new mesh points are produced from old mesh points. In the example below,
32   we begin with an edge (v0, v3). It is extruded for two layers. The original vertex v0 produces two edges, e1 and e2, and three vertices,
33   v0, v2, and v2. Similarly, vertex v3 produces e3, e4, v3, v4, and v5. The original edge produces itself, e5 and e6, as well as face1 and
34   face2. The new mesh points are given the same labels as the original points which produced them. Thus, if v0 had a label value 1, then so
35   would v1, v2, e1 and e2.
36 
37 .vb
38   v2----- e6    -----v5
39   |                  |
40   e2     face2       e4
41   |                  |
42   v1----- e5    -----v4
43   |                  |
44   e1     face1       e3
45   |                  |
46   v0--- original ----v3
47 .ve
48 
49 .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMExtrude()`, `DMPlexTransform`, `DMPlexTransformExtrudeSetThickness()`, `DMPlexTransformExtrudeSetTensor()`
50 @*/
51 PetscErrorCode DMPlexExtrude(DM dm, PetscInt layers, PetscReal thickness, PetscBool tensor, PetscBool symmetric, PetscBool periodic, const PetscReal normal[], const PetscReal thicknesses[], DM *edm)
52 {
53   DMPlexTransform tr;
54   DM              cdm;
55   PetscObject     disc;
56   PetscClassId    id;
57   const char     *prefix;
58   PetscOptions    options;
59   PetscBool       useCeed;
60   PetscBool       cutMarker = PETSC_FALSE;
61 
62   PetscFunctionBegin;
63   PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
64   PetscCall(PetscObjectSetName((PetscObject)tr, "Extrusion Transform"));
65   PetscCall(DMPlexTransformSetDM(tr, dm));
66   PetscCall(DMPlexTransformSetType(tr, DMPLEXEXTRUDE));
67   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
68   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
69   PetscCall(PetscObjectGetOptions((PetscObject)dm, &options));
70   PetscCall(PetscObjectSetOptions((PetscObject)tr, options));
71   PetscCall(DMPlexTransformExtrudeSetLayers(tr, layers));
72   if (thickness > 0.) PetscCall(DMPlexTransformExtrudeSetThickness(tr, thickness));
73   PetscCall(DMPlexTransformExtrudeSetTensor(tr, tensor));
74   PetscCall(DMPlexTransformExtrudeSetSymmetric(tr, symmetric));
75   PetscCall(DMPlexTransformExtrudeSetPeriodic(tr, periodic));
76   if (normal) PetscCall(DMPlexTransformExtrudeSetNormal(tr, normal));
77   if (thicknesses) PetscCall(DMPlexTransformExtrudeSetThicknesses(tr, layers, thicknesses));
78   PetscCall(DMPlexTransformSetFromOptions(tr));
79   PetscCall(PetscObjectSetOptions((PetscObject)tr, NULL));
80   PetscCall(DMPlexTransformSetUp(tr));
81   PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_transform_view"));
82   PetscCall(DMPlexTransformApply(tr, dm, edm));
83   PetscCall(DMPlexGetUseCeed(dm, &useCeed));
84   PetscCall(DMPlexSetUseCeed(*edm, useCeed));
85   PetscCall(DMCopyDisc(dm, *edm));
86   // Handle periodic viewing
87   PetscCall(PetscOptionsGetBool(options, ((PetscObject)dm)->prefix, "-dm_plex_periodic_cut", &cutMarker, NULL));
88   PetscCall(DMPlexTransformExtrudeGetPeriodic(tr, &periodic));
89   if (periodic && cutMarker) {
90     DMLabel  cutLabel;
91     PetscInt dim, pStart, pEnd;
92 
93     PetscCall(DMGetDimension(dm, &dim));
94     PetscCall(DMCreateLabel(*edm, "periodic_cut"));
95     PetscCall(DMGetLabel(*edm, "periodic_cut", &cutLabel));
96     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
97     for (PetscInt p = pStart; p < pEnd; ++p) {
98       DMPolytopeType  ct;
99       DMPolytopeType *rct;
100       PetscInt       *rsize, *rcone, *rornt;
101       PetscInt        Nct;
102 
103       PetscCall(DMPlexGetCellType(dm, p, &ct));
104       PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
105       for (PetscInt n = 0; n < Nct; ++n) {
106         PetscInt pNew, pdim = DMPolytopeTypeGetDim(rct[n]);
107 
108         if (ct == rct[n] || pdim > dim) {
109           PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[n], p, 0, &pNew));
110           PetscCall(DMLabelSetValue(cutLabel, pNew, !pdim ? 1 : 2));
111         }
112       }
113     }
114   }
115   // It is too hard to raise the dimension of a discretization, so just remake it
116   PetscCall(DMGetCoordinateDM(dm, &cdm));
117   PetscCall(DMGetField(cdm, 0, NULL, &disc));
118   PetscCall(PetscObjectGetClassId(disc, &id));
119   if (id == PETSCFE_CLASSID) {
120     PetscSpace sp;
121     PetscInt   deg;
122 
123     PetscCall(PetscFEGetBasisSpace((PetscFE)disc, &sp));
124     PetscCall(PetscSpaceGetDegree(sp, &deg, NULL));
125     PetscCall(DMPlexCreateCoordinateSpace(*edm, deg, PETSC_TRUE, NULL));
126   }
127   PetscCall(DMPlexTransformCreateDiscLabels(tr, *edm));
128   PetscCall(DMPlexTransformDestroy(&tr));
129   if (*edm) {
130     ((DM_Plex *)(*edm)->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
131     ((DM_Plex *)(*edm)->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
132   }
133   PetscFunctionReturn(PETSC_SUCCESS);
134 }
135 
136 PetscErrorCode DMExtrude_Plex(DM dm, PetscInt layers, DM *edm)
137 {
138   PetscFunctionBegin;
139   PetscCall(DMPlexExtrude(dm, layers, PETSC_DETERMINE, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, NULL, NULL, edm));
140   PetscCall(DMSetMatType(*edm, dm->mattype));
141   PetscCall(DMViewFromOptions(*edm, NULL, "-check_extrude"));
142   PetscFunctionReturn(PETSC_SUCCESS);
143 }
144