petsc-3.13.6 2020-09-29
Report Typos and Errors

DMPlexExtrude

Creates a (d+1)-D mesh by extruding a d-D mesh in the normal direction using prismatic cells.

Synopsis

#include "petscdmplex.h"   
#include "petscdmplex.h"   
PetscErrorCode DMPlexExtrude(DM idm, PetscInt layers, PetscReal height, PetscBool ordExt, PetscBool interpolate, DM* dm)
Collective on idm

Input Parameters

idm - The mesh to be extruted
layers - The number of layers
height - The height of the extruded layer
ordExt - If PETSC_TRUE, orders the extruded cells in the height first. Otherwise, orders the cell on the layers first
interpolate - Flag to create intermediate mesh pieces (edges, faces)

Output Parameter

dm -The DM object

Notes: The mesh created has prismatic cells, and the vertex ordering in the cone of the cell is that of the tensor prismatic cells.

See Also

DMPlexCreateWedgeCylinderMesh(), DMPlexCreateWedgeBoxMesh(), DMSetType(), DMCreate()

Level

advanced

Location

src/dm/impls/plex/plexcreate.c
Index of all DMPLEX routines
Table of Contents for all manual pages
Index of all manual pages