DMDACreate3d#
Creates an object that will manage the communication of three-dimensional regular array data that is distributed across one or more MPI processes.
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
Input Parameters#
comm - MPI communicator
bx - type of x ghost nodes the array have. Use one of
DM_BOUNDARY_NONE
,DM_BOUNDARY_GHOSTED
,DM_BOUNDARY_PERIODIC
.by - type of y ghost nodes the array have. Use one of
DM_BOUNDARY_NONE
,DM_BOUNDARY_GHOSTED
,DM_BOUNDARY_PERIODIC
.bz - type of z 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
orDMDA_STENCIL_BOX
)M - global dimension in x direction of the array
N - global dimension in y direction of the array
P - global dimension in z 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)p - corresponding number of processors in z 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
.lz - arrays containing the number of nodes in each cell along the z coordinates, or
NULL
.
Output Parameter#
da - the resulting distributed array object
Options Database Keys#
-dm_view - Calls
DMView()
at the conclusion ofDMDACreate3d()
-da_grid_x
- number of grid points in x direction-da_grid_y
- number of grid points in y direction-da_grid_z
- number of grid points in z direction-da_processors_x
- number of processors in x direction-da_processors_y
- number of processors in y direction-da_processors_z
- number of processors in z direction-da_bd_x
- boundary type in x direction-da_bd_y
- boundary type in y direction-da_bd_z
- boundary type in x 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_z
- refinement ratio in z directio-da_refine
- refine theDMDA
n times before creating it
Notes#
If lx
, ly
, or lz
are 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
.
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.
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()
, DMDACreate2d()
, DMGlobalToLocalBegin()
, DMDAGetRefinementFactor()
,
DMGlobalToLocalEnd()
, DMLocalToGlobalBegin()
, DMLocalToLocalBegin()
, DMLocalToLocalEnd()
, DMDASetRefinementFactor()
,
DMDAGetInfo()
, DMCreateGlobalVector()
, DMCreateLocalVector()
, DMDACreateNaturalVector()
, DMLoad()
, DMDAGetOwnershipRanges()
,
DMStagCreate3d()
, DMBoundaryType
Level#
beginner
Location#
Examples#
src/dm/field/tutorials/ex1.c
src/dm/tutorials/ex11f90.F90
src/dm/tutorials/ex15.c
src/dm/tutorials/ex22.c
src/dm/tutorials/ex14.c
src/dm/tutorials/ex20.c
src/ts/tutorials/ex14.c
src/ksp/ksp/tutorials/ex59.c
src/dm/tutorials/ex13f90.F90
src/dm/tutorials/ex3.c
Index of all DMDA routines
Table of Contents for all manual pages
Index of all manual pages