xref: /petsc/src/dm/impls/plex/plexextrude.c (revision 0619917b5a674bb687c64e7daba2ab22be99af31)
1 #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
2 #include <petscdmplextransform.h>
3 
4 /*@C
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 
60   PetscFunctionBegin;
61   PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
62   PetscCall(DMPlexTransformSetDM(tr, dm));
63   PetscCall(DMPlexTransformSetType(tr, DMPLEXEXTRUDE));
64   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
65   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
66   PetscCall(PetscObjectGetOptions((PetscObject)dm, &options));
67   PetscCall(PetscObjectSetOptions((PetscObject)tr, options));
68   PetscCall(DMPlexTransformExtrudeSetLayers(tr, layers));
69   if (thickness > 0.) PetscCall(DMPlexTransformExtrudeSetThickness(tr, thickness));
70   PetscCall(DMPlexTransformExtrudeSetTensor(tr, tensor));
71   PetscCall(DMPlexTransformExtrudeSetSymmetric(tr, symmetric));
72   PetscCall(DMPlexTransformExtrudeSetPeriodic(tr, periodic));
73   if (normal) PetscCall(DMPlexTransformExtrudeSetNormal(tr, normal));
74   if (thicknesses) PetscCall(DMPlexTransformExtrudeSetThicknesses(tr, layers, thicknesses));
75   PetscCall(DMPlexTransformSetFromOptions(tr));
76   PetscCall(PetscObjectSetOptions((PetscObject)tr, NULL));
77   PetscCall(DMPlexTransformSetUp(tr));
78   PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_transform_view"));
79   PetscCall(DMPlexTransformApply(tr, dm, edm));
80   PetscCall(DMCopyDisc(dm, *edm));
81   // It is too hard to raise the dimension of a discretization, so just remake it
82   PetscCall(DMGetCoordinateDM(dm, &cdm));
83   PetscCall(DMGetField(cdm, 0, NULL, &disc));
84   PetscCall(PetscObjectGetClassId(disc, &id));
85   if (id == PETSCFE_CLASSID) {
86     PetscSpace sp;
87     PetscInt   deg;
88 
89     PetscCall(PetscFEGetBasisSpace((PetscFE)disc, &sp));
90     PetscCall(PetscSpaceGetDegree(sp, &deg, NULL));
91     PetscCall(DMPlexCreateCoordinateSpace(*edm, deg, NULL));
92   }
93   PetscCall(DMPlexTransformCreateDiscLabels(tr, *edm));
94   PetscCall(DMPlexTransformDestroy(&tr));
95   if (*edm) {
96     ((DM_Plex *)(*edm)->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
97     ((DM_Plex *)(*edm)->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
98   }
99   PetscFunctionReturn(PETSC_SUCCESS);
100 }
101 
102 PetscErrorCode DMExtrude_Plex(DM dm, PetscInt layers, DM *edm)
103 {
104   PetscFunctionBegin;
105   PetscCall(DMPlexExtrude(dm, layers, PETSC_DETERMINE, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, NULL, NULL, edm));
106   PetscCall(DMSetMatType((*edm), dm->mattype));
107   PetscCall(DMViewFromOptions(*edm, NULL, "-check_extrude"));
108   PetscFunctionReturn(PETSC_SUCCESS);
109 }
110