DMDA
The DMDA (Distributed Array) module provides functionality for creating and managing structured grids in 1D, 2D, or 3D.
Overview
DMDA is ideal for problems on regular structured grids where:
- The grid is logically rectangular
- Each grid point has the same number of degrees of freedom
- Stencil operations follow a regular pattern (star or box stencils)
Creating a DMDA
# 2D grid example
da = DMDA(
petsclib,
MPI.COMM_WORLD,
(DM_BOUNDARY_NONE, DM_BOUNDARY_NONE), # boundary types
(nx, ny), # global dimensions
1, # degrees of freedom per node
1, # stencil width
DMDA_STENCIL_STAR # stencil type
)Functions
PETSc.DMDA — Method
DMDA(
petsclib::PetscLib
comm::MPI.Comm,
boundary_type::NTuple{D, DMBoundaryType},
global_dim::NTuple{D, Integer},
dof_per_node::Integer,
stencil_width::Integer,
stencil_type;
points_per_proc::Tuple,
processors::Tuple,
setfromoptions = true,
dmsetup = true,
prefix = "",
options...
)Creates a D-dimensional distributed array with the options specified using keyword arguments.
If keyword argument points_per_proc[k] isa Vector{petsclib.PetscInt} then this specifies the points per processor in dimension k.
If keyword argument processors[k] isa Integer then this specifies the number of processors used in dimension k; ignored when D == 1.
If keyword argument setfromoptions == true then setfromoptions! called.
If keyword argument dmsetup == true then setup! is called.
When D == 1 the stencil_type argument is not required and ignored if specified.
External Links
- PETSc Manual:
DMDA/DMDACreate1d
- PETSc Manual:
DMDA/DMDACreate2d
- PETSc Manual:
DMDA/DMDACreate3d
PETSc.localinteriorlinearindex — Method
ind = localinteriorlinearindex(dmda::AbstractPetscDM)Returns the linear indices associated with the degrees of freedom own by this MPI rank embedded in the ghost index space for the dmda
PETSc.ndofs — Method
ndofs(da::AbstractPetscDM)Return the number of dofs in for da
External Links
- PETSc Manual:
DMDA/DMDAGetDof
PETSc.reshapelocalarray — Method
reshapelocalarray(Arr, da::AbstractPetscDM{PetscLib}, ndof = ndofs(da))Returns an array with the same data as Arr but reshaped as an array that can be addressed with global indexing.