petsc-3.7.1 2016-05-15
Report Typos and Errors

DMDACreate2d

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

Synopsis

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

Input Parameters

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

Output Parameter

da -the resulting distributed array object Many br

Options Database Key

-dm_view - Calls DMView() at the conclusion of DMDACreate2d() Many br
-da_grid_x <nx> - number of grid points in x direction, if M < 0 Many br
-da_grid_y <ny> - number of grid points in y direction, if N < 0 Many br
-da_processors_x <nx> - number of processors in x direction Many br
-da_processors_y <ny> - number of processors in y direction Many br
-da_refine_x <rx> - refinement ratio in x direction Many br
-da_refine_y <ry> - refinement ratio in y direction Many br
-da_refine <n> - refine the DMDA n times before creating, if M or N < 0 Many br

Many br

Notes

The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the Many brstandard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes Many brthe standard 9-pt stencil. Many br

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

Keywords

distributed array, create, two-dimensional

See Also

DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(),
DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(), Many brDMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges() Many br

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

Examples

src/dm/examples/tutorials/ex1.c.html
src/dm/examples/tutorials/ex2.c.html
src/dm/examples/tutorials/ex3.c.html
src/dm/examples/tutorials/ex4.c.html
src/dm/examples/tutorials/ex5.c.html
src/dm/examples/tutorials/ex7.c.html
src/dm/examples/tutorials/ex9.c.html
src/dm/examples/tutorials/ex10.c.html
src/dm/examples/tutorials/ex12.c.html
src/dm/examples/tutorials/ex51.c.html
src/dm/examples/tutorials/ex65dm.c.html