petsc-3.5.4 2015-05-23
Report Typos and Errors

DMDACreate3d

Creates an object that will manage the communication of three-dimensional regular array data that is distributed across some processors.

Synopsis

#include "petscdmda.h"    
PetscErrorCode  DMDACreate3d(MPI_Comm comm,DMBoundaryType bx,DMBoundaryType by,DMBoundaryType bz,DMDAStencilType stencil_type,PetscInt M,
               PetscInt N,PetscInt P,PetscInt m,PetscInt n,PetscInt p,PetscInt dof,PetscInt s,const PetscInt lx[],const PetscInt ly[],const PetscInt lz[],DM *da)
Collective on MPI_Comm

Input Parameters

comm - MPI communicator
bx,by,bz - type of ghost nodes the array have. Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.
stencil_type - Type of stencil (DMDA_STENCIL_STAR or DMDA_STENCIL_BOX)
M,N,P - global dimension in each direction of the array (use -M, -N, and or -P to indicate that it may be set to a different value from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>)
m,n,p - corresponding number of processors in each dimension (or PETSC_DECIDE to have calculated)
dof - number of degrees of freedom per node
s - stencil width
lx, ly, lz - arrays containing the number of nodes in each cell along the x, y, and z coordinates, or NULL. If non-null, these must be of length as m,n,p and the corresponding m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of the ly[] must N, sum of the lz[] must be P

Output Parameter

da -the resulting distributed array object

Options Database Key

-dm_view - Calls DMView() at the conclusion of DMDACreate3d()
-da_grid_x <nx> - number of grid points in x direction, if M < 0
-da_grid_y <ny> - number of grid points in y direction, if N < 0
-da_grid_z <nz> - number of grid points in z direction, if P < 0
-da_processors_x <MX> - number of processors in x direction
-da_processors_y <MY> - number of processors in y direction
-da_processors_z <MZ> - number of processors in z direction
-da_refine_x <rx> - refinement ratio in x direction
-da_refine_y <ry> - refinement ratio in y direction
-da_refine_z <rz>- refinement ratio in z directio
-da_refine <n> - refine the DMDA n times before creating it, , if M, N, or P < 0

Notes

The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the standard 7-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes the standard 27-pt stencil.

The array data itself is NOT stored in the DMDA, it is stored in Vec objects; The appropriate vector objects can be obtained with calls to DMCreateGlobalVector() and DMCreateLocalVector() and calls to VecDuplicate() if more are needed.

Keywords

distributed array, create, three-dimensional

See Also

DMDestroy(), DMView(), DMDACreate1d(), DMDACreate2d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(), DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges()

Level:beginner
Location:
src/dm/impls/da/da3.c
Index of all DM routines
Table of Contents for all manual pages
Index of all manual pages

Examples

src/dm/examples/tutorials/ex3.c.html
src/dm/examples/tutorials/ex15.c.html
src/dm/examples/tutorials/ex11f90.F.html
src/ksp/ksp/examples/tutorials/ex34.c.html
src/ksp/ksp/examples/tutorials/ex42.c.html
src/ksp/ksp/examples/tutorials/ex45.c.html
src/ksp/ksp/examples/tutorials/ex22f.F.html
src/snes/examples/tutorials/ex14.c.html
src/snes/examples/tutorials/ex20.c.html
src/snes/examples/tutorials/ex48.c.html
src/snes/examples/tutorials/ex633d_db.c.html