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 131fcf445aSMatthew G. Knepley . periodic - Flag to extrude periodically 149695643eSMatthew G. Knepley . normal - Surface normal vector, or NULL 159695643eSMatthew G. Knepley - thicknesses - Thickness of each layer, or NULL 169695643eSMatthew G. Knepley 179695643eSMatthew G. Knepley Output Parameter: 189695643eSMatthew G. Knepley . edm - The volumetric mesh 199695643eSMatthew G. Knepley 20a1cb98faSBarry Smith Options Database Keys: 219695643eSMatthew G. Knepley + -dm_plex_transform_extrude_thickness <t> - The total thickness of extruded layers 229695643eSMatthew G. Knepley . -dm_plex_transform_extrude_use_tensor <bool> - Use tensor cells when extruding 239695643eSMatthew G. Knepley . -dm_plex_transform_extrude_symmetric <bool> - Extrude layers symmetrically about the surface 241fcf445aSMatthew G. Knepley . -dm_plex_transform_extrude_periodic <bool> - Extrude layers periodically 259695643eSMatthew G. Knepley . -dm_plex_transform_extrude_normal <n0,...,nd> - Specify the extrusion direction 269695643eSMatthew G. Knepley - -dm_plex_transform_extrude_thicknesses <t0,...,tl> - Specify thickness of each layer 279695643eSMatthew G. Knepley 289695643eSMatthew G. Knepley Level: intermediate 299695643eSMatthew G. Knepley 30a1cb98faSBarry Smith Notes: 31da81f932SPierre Jolivet Extrusion is implemented as a `DMPlexTransform`, so that new mesh points are produced from old mesh points. In the example below, 32a1cb98faSBarry 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, 33a1cb98faSBarry 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 34a1cb98faSBarry 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 35a1cb98faSBarry Smith would v1, v2, e1 and e2. 36a1cb98faSBarry Smith 37a1cb98faSBarry Smith .vb 38a1cb98faSBarry Smith v2----- e6 -----v5 39a1cb98faSBarry Smith | | 40a1cb98faSBarry Smith e2 face2 e4 41a1cb98faSBarry Smith | | 42a1cb98faSBarry Smith v1----- e5 -----v4 43a1cb98faSBarry Smith | | 44a1cb98faSBarry Smith e1 face1 e3 45a1cb98faSBarry Smith | | 46a1cb98faSBarry Smith v0--- original ----v3 47a1cb98faSBarry Smith .ve 48a1cb98faSBarry Smith 491cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMExtrude()`, `DMPlexTransform`, `DMPlexTransformExtrudeSetThickness()`, `DMPlexTransformExtrudeSetTensor()` 509695643eSMatthew G. Knepley @*/ 511fcf445aSMatthew G. Knepley PetscErrorCode DMPlexExtrude(DM dm, PetscInt layers, PetscReal thickness, PetscBool tensor, PetscBool symmetric, PetscBool periodic, const PetscReal normal[], const PetscReal thicknesses[], DM *edm) 52d71ae5a4SJacob Faibussowitsch { 53d410b0cfSMatthew G. Knepley DMPlexTransform tr; 541a55e319SMatthew G. Knepley DM cdm; 551a55e319SMatthew G. Knepley PetscObject disc; 561a55e319SMatthew G. Knepley PetscClassId id; 57d410b0cfSMatthew G. Knepley const char *prefix; 58d410b0cfSMatthew G. Knepley PetscOptions options; 59d2b2dc1eSMatthew G. Knepley PetscBool useCeed; 60afdb71dcSMatthew G. Knepley PetscBool cutMarker = PETSC_FALSE; 61d410b0cfSMatthew G. Knepley 62d410b0cfSMatthew G. Knepley PetscFunctionBegin; 639566063dSJacob Faibussowitsch PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr)); 649566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetDM(tr, dm)); 659566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetType(tr, DMPLEXEXTRUDE)); 669566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 679566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix)); 689566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptions((PetscObject)dm, &options)); 699566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptions((PetscObject)tr, options)); 709566063dSJacob Faibussowitsch PetscCall(DMPlexTransformExtrudeSetLayers(tr, layers)); 719566063dSJacob Faibussowitsch if (thickness > 0.) PetscCall(DMPlexTransformExtrudeSetThickness(tr, thickness)); 729566063dSJacob Faibussowitsch PetscCall(DMPlexTransformExtrudeSetTensor(tr, tensor)); 739566063dSJacob Faibussowitsch PetscCall(DMPlexTransformExtrudeSetSymmetric(tr, symmetric)); 7480cd373eSmarkadams4 PetscCall(DMPlexTransformExtrudeSetPeriodic(tr, periodic)); 759566063dSJacob Faibussowitsch if (normal) PetscCall(DMPlexTransformExtrudeSetNormal(tr, normal)); 769566063dSJacob Faibussowitsch if (thicknesses) PetscCall(DMPlexTransformExtrudeSetThicknesses(tr, layers, thicknesses)); 779566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetFromOptions(tr)); 789566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptions((PetscObject)tr, NULL)); 799566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetUp(tr)); 809566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_transform_view")); 819566063dSJacob Faibussowitsch PetscCall(DMPlexTransformApply(tr, dm, edm)); 82d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexGetUseCeed(dm, &useCeed)); 83d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexSetUseCeed(*edm, useCeed)); 849566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, *edm)); 85afdb71dcSMatthew G. Knepley // Handle periodic viewing 86afdb71dcSMatthew G. Knepley PetscCall(PetscOptionsGetBool(options, ((PetscObject)dm)->prefix, "-dm_plex_periodic_cut", &cutMarker, NULL)); 87afdb71dcSMatthew G. Knepley PetscCall(DMPlexTransformExtrudeGetPeriodic(tr, &periodic)); 88afdb71dcSMatthew G. Knepley if (periodic && cutMarker) { 89afdb71dcSMatthew G. Knepley DMLabel cutLabel; 90afdb71dcSMatthew G. Knepley PetscInt dim, pStart, pEnd; 91afdb71dcSMatthew G. Knepley 92afdb71dcSMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 93afdb71dcSMatthew G. Knepley PetscCall(DMCreateLabel(*edm, "periodic_cut")); 94afdb71dcSMatthew G. Knepley PetscCall(DMGetLabel(*edm, "periodic_cut", &cutLabel)); 95afdb71dcSMatthew G. Knepley PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 96afdb71dcSMatthew G. Knepley for (PetscInt p = pStart; p < pEnd; ++p) { 97afdb71dcSMatthew G. Knepley DMPolytopeType ct; 98afdb71dcSMatthew G. Knepley DMPolytopeType *rct; 99afdb71dcSMatthew G. Knepley PetscInt *rsize, *rcone, *rornt; 100afdb71dcSMatthew G. Knepley PetscInt Nct; 101afdb71dcSMatthew G. Knepley 102afdb71dcSMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, p, &ct)); 103afdb71dcSMatthew G. Knepley PetscCall(DMPlexTransformCellTransform(tr, ct, p, NULL, &Nct, &rct, &rsize, &rcone, &rornt)); 104afdb71dcSMatthew G. Knepley for (PetscInt n = 0; n < Nct; ++n) { 105afdb71dcSMatthew G. Knepley PetscInt pNew, pdim = DMPolytopeTypeGetDim(rct[n]); 106afdb71dcSMatthew G. Knepley 107afdb71dcSMatthew G. Knepley if (ct == rct[n] || pdim > dim) { 108afdb71dcSMatthew G. Knepley PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[n], p, 0, &pNew)); 109afdb71dcSMatthew G. Knepley PetscCall(DMLabelSetValue(cutLabel, pNew, !pdim ? 1 : 2)); 110afdb71dcSMatthew G. Knepley } 111afdb71dcSMatthew G. Knepley } 112afdb71dcSMatthew G. Knepley } 113afdb71dcSMatthew G. Knepley } 1141a55e319SMatthew G. Knepley // It is too hard to raise the dimension of a discretization, so just remake it 1159566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 1161a55e319SMatthew G. Knepley PetscCall(DMGetField(cdm, 0, NULL, &disc)); 1171a55e319SMatthew G. Knepley PetscCall(PetscObjectGetClassId(disc, &id)); 1181a55e319SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 1191a55e319SMatthew G. Knepley PetscSpace sp; 1201a55e319SMatthew G. Knepley PetscInt deg; 1211a55e319SMatthew G. Knepley 1221a55e319SMatthew G. Knepley PetscCall(PetscFEGetBasisSpace((PetscFE)disc, &sp)); 1231a55e319SMatthew G. Knepley PetscCall(PetscSpaceGetDegree(sp, °, NULL)); 124e44f6aebSMatthew G. Knepley PetscCall(DMPlexCreateCoordinateSpace(*edm, deg, PETSC_TRUE, NULL)); 1251a55e319SMatthew G. Knepley } 1269566063dSJacob Faibussowitsch PetscCall(DMPlexTransformCreateDiscLabels(tr, *edm)); 1279566063dSJacob Faibussowitsch PetscCall(DMPlexTransformDestroy(&tr)); 128d410b0cfSMatthew G. Knepley if (*edm) { 129d410b0cfSMatthew G. Knepley ((DM_Plex *)(*edm)->data)->printFEM = ((DM_Plex *)dm->data)->printFEM; 130d410b0cfSMatthew G. Knepley ((DM_Plex *)(*edm)->data)->printL2 = ((DM_Plex *)dm->data)->printL2; 131d410b0cfSMatthew G. Knepley } 1323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 133d410b0cfSMatthew G. Knepley } 134d410b0cfSMatthew G. Knepley 135d71ae5a4SJacob Faibussowitsch PetscErrorCode DMExtrude_Plex(DM dm, PetscInt layers, DM *edm) 136d71ae5a4SJacob Faibussowitsch { 137d410b0cfSMatthew G. Knepley PetscFunctionBegin; 1381fcf445aSMatthew G. Knepley PetscCall(DMPlexExtrude(dm, layers, PETSC_DETERMINE, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, NULL, NULL, edm)); 139*f4f49eeaSPierre Jolivet PetscCall(DMSetMatType(*edm, dm->mattype)); 1409566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*edm, NULL, "-check_extrude")); 1413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 142d410b0cfSMatthew G. Knepley } 143