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 109695643eSMatthew G. Knepley . 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 199695643eSMatthew G. Knepley Notes: 209695643eSMatthew G. Knepley Extrusion is implemented as a DMPlexTransform, so that new mesh points are produced from old mesh points. In the exmaple below, 219695643eSMatthew G. Knepley 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, 229695643eSMatthew G. Knepley 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 239695643eSMatthew G. Knepley 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 249695643eSMatthew G. Knepley would v1, v2, e1 and e2. 259695643eSMatthew G. Knepley 269695643eSMatthew G. Knepley $ v2----- e6 -----v5 279695643eSMatthew G. Knepley $ | | 289695643eSMatthew G. Knepley $ e2 face2 e4 299695643eSMatthew G. Knepley $ | | 309695643eSMatthew G. Knepley $ v1----- e5 -----v4 319695643eSMatthew G. Knepley $ | | 329695643eSMatthew G. Knepley $ e1 face1 e3 339695643eSMatthew G. Knepley $ | | 349695643eSMatthew G. Knepley $ v0--- original ----v3 359695643eSMatthew G. Knepley 369695643eSMatthew G. Knepley Options Database: 379695643eSMatthew G. Knepley + -dm_plex_transform_extrude_thickness <t> - The total thickness of extruded layers 389695643eSMatthew G. Knepley . -dm_plex_transform_extrude_use_tensor <bool> - Use tensor cells when extruding 399695643eSMatthew G. Knepley . -dm_plex_transform_extrude_symmetric <bool> - Extrude layers symmetrically about the surface 409695643eSMatthew G. Knepley . -dm_plex_transform_extrude_normal <n0,...,nd> - Specify the extrusion direction 419695643eSMatthew G. Knepley - -dm_plex_transform_extrude_thicknesses <t0,...,tl> - Specify thickness of each layer 429695643eSMatthew G. Knepley 439695643eSMatthew G. Knepley Level: intermediate 449695643eSMatthew G. Knepley 45db781477SPatrick Sanan .seealso: `DMExtrude()`, `DMPlexTransform`, `DMPlexTransformExtrudeSetThickness()`, `DMPlexTransformExtrudeSetTensor()` 469695643eSMatthew G. Knepley @*/ 47d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexExtrude(DM dm, PetscInt layers, PetscReal thickness, PetscBool tensor, PetscBool symmetric, const PetscReal normal[], const PetscReal thicknesses[], DM *edm) 48d71ae5a4SJacob Faibussowitsch { 49d410b0cfSMatthew G. Knepley DMPlexTransform tr; 50*1a55e319SMatthew G. Knepley DM cdm; 51*1a55e319SMatthew G. Knepley PetscObject disc; 52*1a55e319SMatthew G. Knepley PetscClassId id; 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)); 76*1a55e319SMatthew G. Knepley // It is too hard to raise the dimension of a discretization, so just remake it 779566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 78*1a55e319SMatthew G. Knepley PetscCall(DMGetField(cdm, 0, NULL, &disc)); 79*1a55e319SMatthew G. Knepley PetscCall(PetscObjectGetClassId(disc, &id)); 80*1a55e319SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 81*1a55e319SMatthew G. Knepley PetscSpace sp; 82*1a55e319SMatthew G. Knepley PetscInt deg; 83*1a55e319SMatthew G. Knepley 84*1a55e319SMatthew G. Knepley PetscCall(PetscFEGetBasisSpace((PetscFE)disc, &sp)); 85*1a55e319SMatthew G. Knepley PetscCall(PetscSpaceGetDegree(sp, °, NULL)); 86*1a55e319SMatthew G. Knepley PetscCall(DMPlexCreateCoordinateSpace(*edm, deg, NULL)); 87*1a55e319SMatthew G. Knepley } 889566063dSJacob Faibussowitsch PetscCall(DMPlexTransformCreateDiscLabels(tr, *edm)); 899566063dSJacob Faibussowitsch PetscCall(DMPlexTransformDestroy(&tr)); 90d410b0cfSMatthew G. Knepley if (*edm) { 91d410b0cfSMatthew G. Knepley ((DM_Plex *)(*edm)->data)->printFEM = ((DM_Plex *)dm->data)->printFEM; 92d410b0cfSMatthew G. Knepley ((DM_Plex *)(*edm)->data)->printL2 = ((DM_Plex *)dm->data)->printL2; 93d410b0cfSMatthew G. Knepley } 94d410b0cfSMatthew G. Knepley PetscFunctionReturn(0); 95d410b0cfSMatthew G. Knepley } 96d410b0cfSMatthew G. Knepley 97d71ae5a4SJacob Faibussowitsch PetscErrorCode DMExtrude_Plex(DM dm, PetscInt layers, DM *edm) 98d71ae5a4SJacob Faibussowitsch { 99d410b0cfSMatthew G. Knepley PetscFunctionBegin; 1009566063dSJacob Faibussowitsch PetscCall(DMPlexExtrude(dm, layers, PETSC_DETERMINE, PETSC_TRUE, PETSC_FALSE, NULL, NULL, edm)); 1019566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*edm, NULL, "-check_extrude")); 102d410b0cfSMatthew G. Knepley PetscFunctionReturn(0); 103d410b0cfSMatthew G. Knepley } 104