DMStagCreate3d#
Create an object to manage data living on the elements, faces, edges, and vertices of a parallelized regular 3D grid.
Synopsis#
PETSC_EXTERN PetscErrorCode DMStagCreate3d(MPI_Comm comm, DMBoundaryType bndx, DMBoundaryType bndy, DMBoundaryType bndz, PetscInt M, PetscInt N, PetscInt P, PetscInt m, PetscInt n, PetscInt p, PetscInt dof0, PetscInt dof1, PetscInt dof2, PetscInt dof3, DMStagStencilType stencilType, PetscInt stencilWidth, const PetscInt lx[], const PetscInt ly[], const PetscInt lz[], DM *dm)
Collective
Input Parameters#
comm - MPI communicator
bndx - x boundary type,
DM_BOUNDARY_NONE
,DM_BOUNDARY_PERIODIC
, orDM_BOUNDARY_GHOSTED
bndy - y boundary type,
DM_BOUNDARY_NONE
,DM_BOUNDARY_PERIODIC
, orDM_BOUNDARY_GHOSTED
bndz - z boundary type,
DM_BOUNDARY_NONE
,DM_BOUNDARY_PERIODIC
, orDM_BOUNDARY_GHOSTED
M - global number of elements in x direction
N - global number of elements in y direction
P - global number of elements in z direction
m - number of ranks in the x direction (may be
PETSC_DECIDE
)n - number of ranks in the y direction (may be
PETSC_DECIDE
)p - number of ranks in the z direction (may be
PETSC_DECIDE
)dof0 - number of degrees of freedom per vertex/0-cell
dof1 - number of degrees of freedom per edge/1-cell
dof2 - number of degrees of freedom per face/2-cell
dof3 - number of degrees of freedom per element/3-cell
stencilType - ghost/halo region type:
DMSTAG_STENCIL_NONE
,DMSTAG_STENCIL_BOX
, orDMSTAG_STENCIL_STAR
stencilWidth - width, in elements, of halo/ghost region
lx - array of local x element counts, of length equal to
m
, summing toM
, orNULL
ly - arrays of local y element counts, of length equal to
n
, summing toN
, orNULL
lz - arrays of local z element counts, of length equal to
p
, summing toP
, orNULL
Output Parameter#
dm - the new
DMSTAG
object
Options Database Keys#
-dm_view - calls
DMViewFromOptions()
at the conclusion ofDMSetUp()
-stag_grid_x
- number of elements in the x direction-stag_grid_y
- number of elements in the y direction-stag_grid_z
- number of elements in the z direction-stag_ranks_x
- number of ranks in the x direction-stag_ranks_y
- number of ranks in the y direction-stag_ranks_z
- number of ranks in the z direction-stag_ghost_stencil_width - width of ghost region, in elements
-stag_boundary_type x <none,ghosted,periodic> -
DMBoundaryType
value-stag_boundary_type y <none,ghosted,periodic> -
DMBoundaryType
value-stag_boundary_type z <none,ghosted,periodic> -
DMBoundaryType
value
Notes#
You must call DMSetUp()
after this call before using the DM
.
If you wish to use the options database (see the keys above) to change values in the DMSTAG
, you must call
DMSetFromOptions()
after this function but before DMSetUp()
.
See Also#
DMSTAG: Staggered, Structured Grid, DMSTAG
, DMStagCreate1d()
, DMStagCreate2d()
, DMDestroy()
, DMView()
, DMCreateGlobalVector()
, DMCreateLocalVector()
, DMLocalToGlobalBegin()
, DMDACreate3d()
Level#
beginner
Location#
Examples#
src/dm/impls/stag/tutorials/ex4.c
src/dm/impls/stag/tutorials/ex3.c
src/dm/impls/stag/tutorials/ex6.c
Index of all DMStag routines
Table of Contents for all manual pages
Index of all manual pages