DMPlexCreateBoxMesh#
Creates a mesh on the tensor product of unit intervals (box) using simplices or tensor cells (hexahedra).
Synopsis#
#include "petscdmplex.h"
#include "petscdmplextransform.h"
PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate, PetscInt localizationHeight, PetscBool sparseLocalize, DM *dm)
Collective
Input Parameters#
comm - The communicator for the
DM
objectdim - The spatial dimension
simplex -
PETSC_TRUE
for simplices,PETSC_FALSE
for tensor cellsfaces - Number of faces per dimension, or
NULL
for (1,) in 1D and (2, 2) in 2D and (1, 1, 1) in 3Dlower - 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
forDM_BOUNDARY_NONE
interpolate - Flag to create intermediate mesh pieces (edges, faces)
localizationHeight - Flag to localize edges and faces in addition to cells; only significant for periodic meshes
sparseLocalize - Flag to localize coordinates only for cells near the periodic boundary; only significant for periodic meshes
Output Parameter#
dm - The
DM
object
Note#
To customize this mesh using options, use
DMCreate(comm, &dm);
DMSetType(dm, DMPLEX);
DMSetFromOptions(dm);
and use the options in DMSetFromOptions()
.
Here is the numbering returned for 2 faces in each direction for tensor cells:
10---17---11---18----12
| | |
| | |
20 2 22 3 24
| | |
| | |
7---15----8---16----9
| | |
| | |
19 0 21 1 23
| | |
| | |
4---13----5---14----6
and for simplicial cells
14----8---15----9----16
|\ 5 |\ 7 |
| \ | \ |
13 2 14 3 15
| 4 \ | 6 \ |
| \ | \ |
11----6---12----7----13
|\ |\ |
| \ 1 | \ 3 |
10 0 11 1 12
| 0 \ | 2 \ |
| \ | \ |
8----4----9----5----10
See Also#
DMPlex: Unstructured Grids, DM
, DMPLEX
, DMSetFromOptions()
, DMPlexCreateFromFile()
, DMPlexCreateHexCylinderMesh()
, DMSetType()
, DMCreate()
Level#
beginner
Location#
Examples#
src/snes/tutorials/ex75.c
src/dm/tutorials/ex21.c
src/dm/tutorials/ex20.c
src/dm/dt/dualspace/impls/lagrange/tutorials/ex2.c
src/dm/field/tutorials/ex1.c
src/ts/tutorials/ex52.c
Implementations#
DMPlexCreateBoxMesh_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