DMPlexCreateWedgeBoxMesh#

Creates a 3-D mesh tessellating the (x,y) plane and extruding in the third direction using wedge cells.

Synopsis#

#include "petscdmplex.h"   
#include "petscdmplextransform.h"   
PetscErrorCode DMPlexCreateWedgeBoxMesh(MPI_Comm comm, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool orderHeight, PetscBool interpolate, DM *dm)

Collective

Input Parameters#

  • comm - The communicator for the DM object

  • faces - Number of faces per dimension, or NULL for (1, 1, 1)

  • lower - The lower left corner, or NULL for (0, 0, 0)

  • upper - The upper right corner, or NULL for (1, 1, 1)

  • periodicity - The boundary type for the X,Y,Z direction, or NULL for DM_BOUNDARY_NONE

  • orderHeight - 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

See Also#

DMPlex: Unstructured Grids, DM, DMPLEX, DMPlexCreateHexCylinderMesh(), DMPlexCreateWedgeCylinderMesh(), DMExtrude(), DMPlexCreateBoxMesh(), DMSetType(), DMCreate()

Level#

beginner

Location#

src/dm/impls/plex/plexcreate.c

Implementations#

DMPlexCreateWedgeBoxMesh_Internal() in src/dm/impls/plex/plexcreate.c


Index of all DMPlex routines
Table of Contents for all manual pages
Index of all manual pages