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 10a1cb98faSBarry 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 19a1cb98faSBarry 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 28a1cb98faSBarry Smith Notes: 29da81f932SPierre Jolivet Extrusion is implemented as a `DMPlexTransform`, so that new mesh points are produced from old mesh points. In the example below, 30a1cb98faSBarry 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, 31a1cb98faSBarry 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 32a1cb98faSBarry 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 33a1cb98faSBarry Smith would v1, v2, e1 and e2. 34a1cb98faSBarry Smith 35a1cb98faSBarry Smith .vb 36a1cb98faSBarry Smith v2----- e6 -----v5 37a1cb98faSBarry Smith | | 38a1cb98faSBarry Smith e2 face2 e4 39a1cb98faSBarry Smith | | 40a1cb98faSBarry Smith v1----- e5 -----v4 41a1cb98faSBarry Smith | | 42a1cb98faSBarry Smith e1 face1 e3 43a1cb98faSBarry Smith | | 44a1cb98faSBarry Smith v0--- original ----v3 45a1cb98faSBarry Smith .ve 46a1cb98faSBarry Smith 47*1cc06b55SBarry Smith .seealso: [](ch_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; 521a55e319SMatthew G. Knepley DM cdm; 531a55e319SMatthew G. Knepley PetscObject disc; 541a55e319SMatthew G. Knepley PetscClassId id; 55d410b0cfSMatthew G. Knepley const char *prefix; 56d410b0cfSMatthew G. Knepley PetscOptions options; 57d410b0cfSMatthew G. Knepley 58d410b0cfSMatthew G. Knepley PetscFunctionBegin; 599566063dSJacob Faibussowitsch PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr)); 609566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetDM(tr, dm)); 619566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetType(tr, DMPLEXEXTRUDE)); 629566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 639566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix)); 649566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptions((PetscObject)dm, &options)); 659566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptions((PetscObject)tr, options)); 669566063dSJacob Faibussowitsch PetscCall(DMPlexTransformExtrudeSetLayers(tr, layers)); 679566063dSJacob Faibussowitsch if (thickness > 0.) PetscCall(DMPlexTransformExtrudeSetThickness(tr, thickness)); 689566063dSJacob Faibussowitsch PetscCall(DMPlexTransformExtrudeSetTensor(tr, tensor)); 699566063dSJacob Faibussowitsch PetscCall(DMPlexTransformExtrudeSetSymmetric(tr, symmetric)); 709566063dSJacob Faibussowitsch if (normal) PetscCall(DMPlexTransformExtrudeSetNormal(tr, normal)); 719566063dSJacob Faibussowitsch if (thicknesses) PetscCall(DMPlexTransformExtrudeSetThicknesses(tr, layers, thicknesses)); 729566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetFromOptions(tr)); 739566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptions((PetscObject)tr, NULL)); 749566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetUp(tr)); 759566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_transform_view")); 769566063dSJacob Faibussowitsch PetscCall(DMPlexTransformApply(tr, dm, edm)); 779566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, *edm)); 781a55e319SMatthew G. Knepley // It is too hard to raise the dimension of a discretization, so just remake it 799566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 801a55e319SMatthew G. Knepley PetscCall(DMGetField(cdm, 0, NULL, &disc)); 811a55e319SMatthew G. Knepley PetscCall(PetscObjectGetClassId(disc, &id)); 821a55e319SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 831a55e319SMatthew G. Knepley PetscSpace sp; 841a55e319SMatthew G. Knepley PetscInt deg; 851a55e319SMatthew G. Knepley 861a55e319SMatthew G. Knepley PetscCall(PetscFEGetBasisSpace((PetscFE)disc, &sp)); 871a55e319SMatthew G. Knepley PetscCall(PetscSpaceGetDegree(sp, °, NULL)); 881a55e319SMatthew G. Knepley PetscCall(DMPlexCreateCoordinateSpace(*edm, deg, NULL)); 891a55e319SMatthew G. Knepley } 909566063dSJacob Faibussowitsch PetscCall(DMPlexTransformCreateDiscLabels(tr, *edm)); 919566063dSJacob Faibussowitsch PetscCall(DMPlexTransformDestroy(&tr)); 92d410b0cfSMatthew G. Knepley if (*edm) { 93d410b0cfSMatthew G. Knepley ((DM_Plex *)(*edm)->data)->printFEM = ((DM_Plex *)dm->data)->printFEM; 94d410b0cfSMatthew G. Knepley ((DM_Plex *)(*edm)->data)->printL2 = ((DM_Plex *)dm->data)->printL2; 95d410b0cfSMatthew G. Knepley } 963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 97d410b0cfSMatthew G. Knepley } 98d410b0cfSMatthew G. Knepley 99d71ae5a4SJacob Faibussowitsch PetscErrorCode DMExtrude_Plex(DM dm, PetscInt layers, DM *edm) 100d71ae5a4SJacob Faibussowitsch { 101d410b0cfSMatthew G. Knepley PetscFunctionBegin; 1029566063dSJacob Faibussowitsch PetscCall(DMPlexExtrude(dm, layers, PETSC_DETERMINE, PETSC_TRUE, PETSC_FALSE, NULL, NULL, edm)); 1039566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*edm, NULL, "-check_extrude")); 1043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 105d410b0cfSMatthew G. Knepley } 106