xref: /petsc/src/dm/impls/plex/plexextrude.c (revision a1cb98fac0cdf0eb4d3e8a0c8b58f3fe8f800bc6) !
1d410b0cfSMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I      "petscdmplex.h"   I*/
2d410b0cfSMatthew G. Knepley #include <petscdmplextransform.h>
3d410b0cfSMatthew G. Knepley 
49695643eSMatthew G. Knepley /*@C
59695643eSMatthew G. Knepley   DMPlexExtrude - Extrude a volumetric mesh from the input surface mesh
69695643eSMatthew G. Knepley 
79695643eSMatthew G. Knepley   Input Parameters:
89695643eSMatthew G. Knepley + dm          - The surface mesh
99695643eSMatthew G. Knepley . layers      - The number of extruded layers
10*a1cb98faSBarry Smith . thickness   - The total thickness of the extruded layers, or `PETSC_DETERMINE`
119695643eSMatthew G. Knepley . tensor      - Flag to create tensor produt cells
129695643eSMatthew G. Knepley . symmetric   - Flag to extrude symmetrically about the surface
139695643eSMatthew G. Knepley . normal      - Surface normal vector, or NULL
149695643eSMatthew G. Knepley - thicknesses - Thickness of each layer, or NULL
159695643eSMatthew G. Knepley 
169695643eSMatthew G. Knepley   Output Parameter:
179695643eSMatthew G. Knepley . edm - The volumetric mesh
189695643eSMatthew G. Knepley 
19*a1cb98faSBarry Smith   Options Database Keys:
209695643eSMatthew G. Knepley + -dm_plex_transform_extrude_thickness <t>           - The total thickness of extruded layers
219695643eSMatthew G. Knepley . -dm_plex_transform_extrude_use_tensor <bool>       - Use tensor cells when extruding
229695643eSMatthew G. Knepley . -dm_plex_transform_extrude_symmetric <bool>        - Extrude layers symmetrically about the surface
239695643eSMatthew G. Knepley . -dm_plex_transform_extrude_normal <n0,...,nd>      - Specify the extrusion direction
249695643eSMatthew G. Knepley - -dm_plex_transform_extrude_thicknesses <t0,...,tl> - Specify thickness of each layer
259695643eSMatthew G. Knepley 
269695643eSMatthew G. Knepley   Level: intermediate
279695643eSMatthew G. Knepley 
28*a1cb98faSBarry Smith   Notes:
29*a1cb98faSBarry Smith   Extrusion is implemented as a `DMPlexTransform`, so that new mesh points are produced from old mesh points. In the exmaple below,
30*a1cb98faSBarry Smith 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,
31*a1cb98faSBarry Smith 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
32*a1cb98faSBarry Smith 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
33*a1cb98faSBarry Smith would v1, v2, e1 and e2.
34*a1cb98faSBarry Smith 
35*a1cb98faSBarry Smith .vb
36*a1cb98faSBarry Smith   v2----- e6    -----v5
37*a1cb98faSBarry Smith   |                  |
38*a1cb98faSBarry Smith   e2     face2       e4
39*a1cb98faSBarry Smith   |                  |
40*a1cb98faSBarry Smith   v1----- e5    -----v4
41*a1cb98faSBarry Smith   |                  |
42*a1cb98faSBarry Smith   e1     face1       e3
43*a1cb98faSBarry Smith   |                  |
44*a1cb98faSBarry Smith   v0--- original ----v3
45*a1cb98faSBarry Smith .ve
46*a1cb98faSBarry Smith 
47*a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `DMPLEX`, `DMExtrude()`, `DMPlexTransform`, `DMPlexTransformExtrudeSetThickness()`, `DMPlexTransformExtrudeSetTensor()`
489695643eSMatthew G. Knepley @*/
49d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexExtrude(DM dm, PetscInt layers, PetscReal thickness, PetscBool tensor, PetscBool symmetric, const PetscReal normal[], const PetscReal thicknesses[], DM *edm)
50d71ae5a4SJacob Faibussowitsch {
51d410b0cfSMatthew G. Knepley   DMPlexTransform tr;
52d410b0cfSMatthew G. Knepley   DM              cdm, ecdm;
53d410b0cfSMatthew G. Knepley   const char     *prefix;
54d410b0cfSMatthew G. Knepley   PetscOptions    options;
55d410b0cfSMatthew G. Knepley 
56d410b0cfSMatthew G. Knepley   PetscFunctionBegin;
579566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
589566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformSetDM(tr, dm));
599566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformSetType(tr, DMPLEXEXTRUDE));
609566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
619566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
629566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptions((PetscObject)dm, &options));
639566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptions((PetscObject)tr, options));
649566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformExtrudeSetLayers(tr, layers));
659566063dSJacob Faibussowitsch   if (thickness > 0.) PetscCall(DMPlexTransformExtrudeSetThickness(tr, thickness));
669566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformExtrudeSetTensor(tr, tensor));
679566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformExtrudeSetSymmetric(tr, symmetric));
689566063dSJacob Faibussowitsch   if (normal) PetscCall(DMPlexTransformExtrudeSetNormal(tr, normal));
699566063dSJacob Faibussowitsch   if (thicknesses) PetscCall(DMPlexTransformExtrudeSetThicknesses(tr, layers, thicknesses));
709566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformSetFromOptions(tr));
719566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptions((PetscObject)tr, NULL));
729566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformSetUp(tr));
739566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_transform_view"));
749566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformApply(tr, dm, edm));
759566063dSJacob Faibussowitsch   PetscCall(DMCopyDisc(dm, *edm));
769566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
779566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(*edm, &ecdm));
789566063dSJacob Faibussowitsch   PetscCall(DMCopyDisc(cdm, ecdm));
799566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformCreateDiscLabels(tr, *edm));
809566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformDestroy(&tr));
81d410b0cfSMatthew G. Knepley   if (*edm) {
82d410b0cfSMatthew G. Knepley     ((DM_Plex *)(*edm)->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
83d410b0cfSMatthew G. Knepley     ((DM_Plex *)(*edm)->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
84d410b0cfSMatthew G. Knepley   }
85d410b0cfSMatthew G. Knepley   PetscFunctionReturn(0);
86d410b0cfSMatthew G. Knepley }
87d410b0cfSMatthew G. Knepley 
88d71ae5a4SJacob Faibussowitsch PetscErrorCode DMExtrude_Plex(DM dm, PetscInt layers, DM *edm)
89d71ae5a4SJacob Faibussowitsch {
90d410b0cfSMatthew G. Knepley   PetscFunctionBegin;
919566063dSJacob Faibussowitsch   PetscCall(DMPlexExtrude(dm, layers, PETSC_DETERMINE, PETSC_TRUE, PETSC_FALSE, NULL, NULL, edm));
929566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*edm, NULL, "-check_extrude"));
93d410b0cfSMatthew G. Knepley   PetscFunctionReturn(0);
94d410b0cfSMatthew G. Knepley }
95