#include "petscdmplex.h" #include "petscdmplex.h" PetscErrorCode DMPlexCreateBoxMesh(MPI_Comm comm, PetscInt dim, PetscBool simplex, const PetscInt faces[], const PetscReal lower[], const PetscReal upper[], const DMBoundaryType periodicity[], PetscBool interpolate, DM *dm)Collective on MPI_Comm
comm | - The communicator for the DM object | |
dim | - The spatial dimension | |
simplex | - PETSC_TRUE for simplices, PETSC_FALSE for tensor cells | |
faces | - Number of faces per dimension, or NULL for (1,) in 1D and (2, 2) in 2D and (1, 1, 1) in 3D | |
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 | |
interpolate | - Flag to create intermediate mesh pieces (edges, faces) |
-dm_plex_box_lower <x,y,z> | - Specify lower-left-bottom coordinates for the box | |
-dm_plex_box_upper <x,y,z> | - Specify upper-right-top coordinates for the box |
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