DMDACreate2d#

Creates an object that will manage the communication of two-dimensional regular array data that is distributed across one or more MPI processes.

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

Input Parameters#

  • comm - MPI communicator

  • bx - type of ghost nodes the x array have. Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.

  • by - type of ghost nodes the y array have. Use one of DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_PERIODIC.

  • stencil_type - stencil type. Use either DMDA_STENCIL_BOX or DMDA_STENCIL_STAR.

  • M - global dimension in x direction of the array

  • N - global dimension in y direction of the array

  • m - corresponding number of processors in x dimension (or PETSC_DECIDE to have calculated)

  • n - corresponding number of processors in y dimension (or PETSC_DECIDE to have calculated)

  • dof - number of degrees of freedom per node

  • s - stencil width

  • lx - arrays containing the number of nodes in each cell along the x coordinates, or NULL.

  • ly - arrays containing the number of nodes in each cell along the y coordinates, or NULL.

Output Parameter#

  • da - the resulting distributed array object

Options Database Keys#

  • -dm_view - Calls DMView() at the conclusion of DMDACreate2d()

  • -da_grid_x - number of grid points in x direction

  • -da_grid_y - number of grid points in y direction

  • -da_processors_x - number of processors in x direction

  • -da_processors_y - number of processors in y direction

  • -da_bd_x - boundary type in x direction

  • -da_bd_y - boundary type in y direction

  • -da_bd_all - boundary type in all directions

  • -da_refine_x - refinement ratio in x direction

  • -da_refine_y - refinement ratio in y direction

  • -da_refine - refine the DMDA n times before creating

Notes#

If lx or ly are non-null, these must be of length as m and n, and the corresponding m and n cannot be PETSC_DECIDE. The sum of the lx entries must be M, and the sum of the ly entries must be N.

The stencil type DMDA_STENCIL_STAR with width 1 corresponds to the standard 5-pt stencil, while DMDA_STENCIL_BOX with width 1 denotes the standard 9-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.

You must call DMSetUp() after this call before using this DM.

To use the options database to change values in the DMDA call DMSetFromOptions() after this call but before DMSetUp().

See Also#

DMDA - Creating vectors for structured grids, DM, DMDA, DMDestroy(), DMView(), DMDACreate1d(), DMDACreate3d(), DMGlobalToLocalBegin(), DMDAGetRefinementFactor(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDASetRefinementFactor(), DMDAGetInfo(), DMCreateGlobalVector(), DMCreateLocalVector(), DMDACreateNaturalVector(), DMLoad(), DMDAGetOwnershipRanges(), DMStagCreate2d(), DMBoundaryType

Level#

beginner

Location#

src/dm/impls/da/da2.c

Examples#

src/dm/field/tutorials/ex1.c
src/tao/unconstrained/tutorials/eptorsion2.c
src/tao/bound/tutorials/plate2.c
src/dm/tutorials/ex21.c
src/dm/tutorials/ex5.c
src/tao/bound/tutorials/jbearing2.c
src/dm/tutorials/ex1.c
src/dm/tutorials/ex9.c
src/tao/unconstrained/tutorials/minsurf2.c
src/dm/tutorials/ex3.c


Index of all DMDA routines
Table of Contents for all manual pages
Index of all manual pages