DMStag

The DMStag (Staggered Grid DM) module provides data management for staggered grids, commonly used in finite difference/volume methods for fluid dynamics and similar applications.

Overview

DMStag is designed for problems where:

  • Variables live at different grid locations (vertices, edges, faces, cell centers)
  • Staggered grids provide better stability for incompressible flow
  • Multiple degrees of freedom per grid location are needed

Staggered Grid Layout

In a staggered grid, different physical quantities are stored at different locations:

1D:

  • Vertices: Scalar quantities (pressure, temperature)
  • Elements: Flux quantities

2D:

  • Vertices: Corner values
  • Edges: Face-normal velocities (u on vertical edges, v on horizontal edges)
  • Elements: Cell-centered values (pressure)

3D:

  • Vertices: Corner values
  • Edges: Edge-centered values
  • Faces: Face-normal quantities
  • Elements: Cell-centered values

Creating a DMStag

# 2D staggered grid
dm = DMStag(
    petsclib,
    MPI.COMM_WORLD,
    (DM_BOUNDARY_NONE, DM_BOUNDARY_NONE),  # boundary types
    (nx, ny),                               # global dimensions
    (dof_vertex, dof_edge, dof_element),   # DOF at each location
    1,                                      # stencil width
    DMSTAG_STENCIL_BOX                     # stencil type
)

# 3D staggered grid
dm = DMStag(
    petsclib,
    MPI.COMM_WORLD,
    (DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE),
    (nx, ny, nz),
    (dof_vertex, dof_edge, dof_face, dof_element),
    1,
    DMSTAG_STENCIL_BOX
)

Accessing Data

Grid Corners and Sizes

# Get local grid extent (without ghost points)
corners = getcorners_dmstag(dm)
# Returns (lower=CartesianIndex, upper=CartesianIndex, size=Tuple)

# Get local grid extent (with ghost points)  
ghost_corners = getghostcorners_dmstag(dm)

Working with Vectors

# Create global and local vectors
global_vec = DMGlobalVec(dm)
local_vec = DMLocalVec(dm)

# Transfer data between global and local
dm_global_to_local!(dm, global_vec, INSERT_VALUES, local_vec)
dm_local_to_global!(dm, local_vec, ADD_VALUES, global_vec)

Getting Location Indices

# Get indices for accessing specific DOF locations
indices = DMStagGetIndices(dm)
# Use indices to access vertex, edge, face, or element DOFs

Setting Coordinates

# Set uniform coordinates
setuniformcoordinates_stag!(dm, xmin, xmax)           # 1D
setuniformcoordinates_stag!(dm, xmin, xmax, ymin, ymax)  # 2D
setuniformcoordinates_stag!(dm, xmin, xmax, ymin, ymax, zmin, zmax)  # 3D

# Get local coordinate array
coords = getlocalcoordinatearray(dm)

Stencil Types

  • DMSTAG_STENCIL_BOX - Full box stencil (includes diagonals)
  • DMSTAG_STENCIL_STAR - Star stencil (axis-aligned neighbors only)

Example: 2D Stokes Flow Setup

# Create staggered grid for Stokes: velocity on edges, pressure in cells
dm = DMStag(
    petsclib,
    MPI.COMM_WORLD,
    (DM_BOUNDARY_NONE, DM_BOUNDARY_NONE),
    (64, 64),      # 64x64 grid
    (0, 1, 1),     # 0 DOF at vertices, 1 at edges (velocity), 1 in elements (pressure)
    1,
    DMSTAG_STENCIL_BOX
)

setuniformcoordinates_stag!(dm, 0.0, 1.0, 0.0, 1.0)

# Create vectors and matrix
x = DMGlobalVec(dm)
b = DMGlobalVec(dm)
A = DMCreateMatrix(dm)

Functions

PETSc.DMStagMethod
da = DMStag(
    petsclib::PetscLib
    comm::MPI.Comm,
    boundary_type::NTuple{D, DMBoundaryType},
    global_dim::NTuple{D, Integer},
    dof_per_node::NTuple{1 + D, Integer},
    stencil_width::Integer,
    stencil_type;
    points_per_proc::Tuple,
    processors::Tuple,
    setfromoptions = true,
    dmsetup = true,
    prefix = "",
    options...
)

Creates a D-dimensional distributed staggered array with the options specified using keyword arguments.

The Tuple dof_per_node specifies how many degrees of freedom are at all the staggerings in the order:

  • 1D: (vertex, element)
  • 2D: (vertex, edge, element)
  • 3D: (vertex, edge, face, element)

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

source
PETSc.DMStagGetIndicesFunction
DMStagGetIndices(dm::DMStag)

Return indices for the central/vertex nodes of a local array built from the input dm. This takes ghost points into account and provides index ranges for accessing staggered data.

Returns

A NamedTuple with:

  • center: Tuple of ranges (x, y, z) for cell-centered indices
  • vertex: Tuple of ranges (x, y, z) for vertex indices

Note

In Julia, array indices start at 1, whereas PETSc uses 0-based indexing with possibly negative ghost indices. This function handles the conversion automatically.

source
PETSc.getcorners_dmstagMethod
corners = getcorners_dmstag(dm::AbstractPetscDM)

Returns a NamedTuple with the global indices (excluding ghost points) of the lower and upper corners as well as the size. Also included is nextra of the number of extra partial elements in each direction.

External Links

source
PETSc.getghostcorners_dmstagMethod
corners = getghostcorners_dmstag(dm::AbstractPetscDM)

Returns a NamedTuple with the global indices (including ghost points) of the lower and upper corners as well as the size. Also included is nextra of the number of extra partial elements in each direction.

External Links

source
PETSc.setuniformcoordinates_stag!Method
setuniformcoordinates_stag!(
    dm::AbstractDMStag,
    xyzmin::Union{NTuple{1,Int},NTuple{2,Int},NTuple{3,Int}},
    xyzmax::Union{NTuple{1,Int},NTuple{2,Int},NTuple{3,Int}},
)

Sets uniform coordinates for the DMStag dm in the range specified by xyzmin and xyzmax.

source