petsc4py.PETSc.DMStag#

class petsc4py.PETSc.DMStag#

Bases: DM

A DM object representing a “staggered grid” or a structured cell complex.

Enumerations

StencilLocation

Stencil location types.

StencilType

Stencil types.

Methods Summary

VecSplitToDMDA(vec, loc, c)

Return DMDA, Vec from a subgrid of a DMStag, its Vec.

create(dim[, dofs, sizes, boundary_types, ...])

Create a DMDA object.

createCompatibleDMStag(dofs)

Create a compatible DMStag with different DOFs/stratum.

get1dCoordinatecArrays()

Not implemented.

getBoundaryTypes()

Return boundary types in each dimension.

getCorners()

Return starting element index, width and number of partial elements.

getDim()

Return the number of dimensions.

getDof()

Get number of DOFs associated with each stratum of the grid.

getEntriesPerElement()

Return the number of entries per element in the local representation.

getGhostCorners()

Return starting element index and width of local region.

getGlobalSizes()

Return global element counts in each dimension.

getIsFirstRank()

Return whether this process is first in each dimension in the process grid.

getIsLastRank()

Return whether this process is last in each dimension in the process grid.

getLocalSizes()

Return local elementwise sizes in each dimension.

getLocationDof(loc)

Return number of DOFs associated with a given point on the grid.

getLocationSlot(loc, c)

Return index to use in accessing raw local arrays.

getOwnershipRanges()

Return elements per process in each dimension.

getProcSizes()

Return number of processes in each dimension.

getProductCoordinateLocationSlot(loc)

Return slot for use with local product coordinate arrays.

getStencilType()

Return elementwise ghost/halo stencil type.

getStencilWidth()

Return elementwise stencil width.

getVecArray(vec)

Not implemented.

migrateVec(vec, dmTo, vecTo)

Transfer a vector between two DMStag objects.

setBoundaryTypes(boundary_types)

Set the boundary types.

setCoordinateDMType(dmtype)

Set the type to store coordinates.

setDof(dofs)

Set DOFs/stratum.

setGlobalSizes(sizes)

Set global element counts in each dimension.

setOwnershipRanges(ranges)

Set elements per process in each dimension.

setProcSizes(sizes)

Set the number of processes in each dimension in the global process grid.

setStencilType(stenciltype)

Set elementwise ghost/halo stencil type.

setStencilWidth(swidth)

Set elementwise stencil width.

setUniformCoordinates([xmin, xmax, ymin, ...])

Set the coordinates to be a uniform grid..

setUniformCoordinatesExplicit([xmin, xmax, ...])

Set coordinates to be a uniform grid, storing all values.

setUniformCoordinatesProduct([xmin, xmax, ...])

Create uniform coordinates, as a product of 1D arrays.

Attributes Summary

boundary_types

Boundary types in each dimension.

corners

The lower left corner and size of local region in each dimension.

dim

The dimension.

dofs

The number of DOFs associated with each stratum of the grid.

entries_per_element

The number of entries per element in the local representation.

ghost_corners

The lower left corner and size of local region in each dimension.

global_sizes

Global element counts in each dimension.

local_sizes

Local elementwise sizes in each dimension.

proc_sizes

The number of processes in each dimension in the global decomposition.

stencil_type

Stencil type.

stencil_width

Elementwise stencil width.

Methods Documentation

VecSplitToDMDA(vec, loc, c)#

Return DMDA, Vec from a subgrid of a DMStag, its Vec.

Collective.

If a c value of -k is provided, the first k DOFs for that position are extracted, padding with zero values if needed. If a non-negative value is provided, a single DOF is extracted.

Parameters:
  • vec (Vec) – The Vec object.

  • loc (StencilLocation) – Which subgrid to extract.

  • c (int) – Which component to extract.

Return type:

tuple[DMDA, Vec]

Source code at petsc4py/PETSc/DMStag.pyx:790

create(dim, dofs=None, sizes=None, boundary_types=None, stencil_type=None, stencil_width=None, proc_sizes=None, ownership_ranges=None, comm=None, setUp=False)#

Create a DMDA object.

Collective.

Creates an object to manage data living on the elements and vertices / the elements, faces, and vertices / the elements, faces, edges, and vertices of a parallelized regular 1D / 2D / 3D grid.

Parameters:
  • dim (int) – The number of dimensions.

  • dofs (tuple[int, ...] | None) – The number of degrees of freedom per vertex, element (1D); vertex, face, element (2D); or vertex, edge, face, element (3D).

  • sizes (tuple[int, ...] | None) – The number of elements in each dimension.

  • boundary_types (tuple[DM.BoundaryType | int | str | bool, ...] | None) – The boundary types.

  • stencil_type (StencilType | None) – The ghost/halo stencil type.

  • stencil_width (int | None) – The width of the ghost/halo region.

  • proc_sizes (tuple[int, ...] | None) – The number of processes in x, y, z dimensions.

  • ownership_ranges (tuple[Sequence[int], ...] | None) – Local x, y, z element counts, of length equal to proc_sizes, summing to sizes.

  • comm (Comm | None) – MPI communicator, defaults to Sys.getDefaultComm.

  • setUp (bool | None) – Whether to call the setup routine after creating the object.

Return type:

Self

Source code at petsc4py/PETSc/DMStag.pyx:50

createCompatibleDMStag(dofs)#

Create a compatible DMStag with different DOFs/stratum.

Collective.

Parameters:

dofs (tuple[int, ...]) – The number of DOFs on the strata in the new DMStag.

Return type:

DM

Source code at petsc4py/PETSc/DMStag.pyx:766

get1dCoordinatecArrays()#

Not implemented.

Source code at petsc4py/PETSc/DMStag.pyx:832

Return type:

None

getBoundaryTypes()#

Return boundary types in each dimension.

Not collective.

Source code at petsc4py/PETSc/DMStag.pyx:478

Return type:

tuple[str, …]

getCorners()#

Return starting element index, width and number of partial elements.

Not collective.

The returned value is calculated excluding ghost points.

The number of extra partial elements is either 1 or 0. The value is 1 on right, top, and front non-periodic domain (“physical”) boundaries, in the x, y, and z dimensions respectively, and otherwise 0.

Source code at petsc4py/PETSc/DMStag.pyx:363

Return type:

tuple[tuple[int, …], tuple[int, …], tuple[int, …]]

getDim()#

Return the number of dimensions.

Not collective.

Source code at petsc4py/PETSc/DMStag.pyx:310

Return type:

int

getDof()#

Get number of DOFs associated with each stratum of the grid.

Not collective.

See also

DMStagGetDOF

Source code at petsc4py/PETSc/DMStag.pyx:348

Return type:

tuple[int, …]

getEntriesPerElement()#

Return the number of entries per element in the local representation.

Not collective.

This is the natural block size for most local operations.

Source code at petsc4py/PETSc/DMStag.pyx:318

Return type:

int

getGhostCorners()#

Return starting element index and width of local region.

Not collective.

Source code at petsc4py/PETSc/DMStag.pyx:385

Return type:

tuple[tuple[int, …], tuple[int, …]]

getGlobalSizes()#

Return global element counts in each dimension.

Not collective.

Source code at petsc4py/PETSc/DMStag.pyx:417

Return type:

tuple[int, …]

getIsFirstRank()#

Return whether this process is first in each dimension in the process grid.

Not collective.

Source code at petsc4py/PETSc/DMStag.pyx:496

Return type:

tuple[int, …]

getIsLastRank()#

Return whether this process is last in each dimension in the process grid.

Not collective.

Source code at petsc4py/PETSc/DMStag.pyx:512

Return type:

tuple[int, …]

getLocalSizes()#

Return local elementwise sizes in each dimension.

Not collective.

The returned value is calculated excluding ghost points.

Source code at petsc4py/PETSc/DMStag.pyx:400

Return type:

tuple[int, …]

getLocationDof(loc)#

Return number of DOFs associated with a given point on the grid.

Not collective.

Parameters:

loc (StencilLocation) – The grid point.

Return type:

int

Source code at petsc4py/PETSc/DMStag.pyx:721

getLocationSlot(loc, c)#

Return index to use in accessing raw local arrays.

Not collective.

Parameters:
Return type:

int

Source code at petsc4py/PETSc/DMStag.pyx:678

getOwnershipRanges()#

Return elements per process in each dimension.

Not collective.

Source code at petsc4py/PETSc/DMStag.pyx:461

Return type:

tuple[Sequence[int], …]

getProcSizes()#

Return number of processes in each dimension.

Not collective.

Source code at petsc4py/PETSc/DMStag.pyx:432

Return type:

tuple[int, …]

getProductCoordinateLocationSlot(loc)#

Return slot for use with local product coordinate arrays.

Not collective.

Parameters:

loc (StencilLocation) – The grid location.

Return type:

None

Source code at petsc4py/PETSc/DMStag.pyx:701

getStencilType()#

Return elementwise ghost/halo stencil type.

Not collective.

Source code at petsc4py/PETSc/DMStag.pyx:447

Return type:

str

getStencilWidth()#

Return elementwise stencil width.

Not collective.

Source code at petsc4py/PETSc/DMStag.pyx:334

Return type:

int

getVecArray(vec)#

Not implemented.

Source code at petsc4py/PETSc/DMStag.pyx:828

Parameters:

vec (Vec)

Return type:

None

migrateVec(vec, dmTo, vecTo)#

Transfer a vector between two DMStag objects.

Collective.

Currently only implemented to migrate global vectors to global vectors.

Parameters:
  • vec (Vec) – The source vector.

  • dmTo (DM) – The compatible destination object.

  • vecTo (Vec) – The destination vector.

Return type:

None

See also

DMStagMigrateVec

Source code at petsc4py/PETSc/DMStag.pyx:743

setBoundaryTypes(boundary_types)#

Set the boundary types.

Logically collective.

Parameters:

boundary_types (tuple[BoundaryType | int | str | bool, ...]) – Boundary types for one/two/three dimensions.

Return type:

None

Source code at petsc4py/PETSc/DMStag.pyx:201

setCoordinateDMType(dmtype)#

Set the type to store coordinates.

Logically collective.

Parameters:

dmtype (Type) – The type to store coordinates.

Return type:

None

Source code at petsc4py/PETSc/DMStag.pyx:657

setDof(dofs)#

Set DOFs/stratum.

Logically collective.

Parameters:

dofs (tuple[int, ...]) – The number of points per 0-cell (vertex/node), 1-cell (element in 1D, edge in 2D and 3D), 2-cell (element in 2D, face in 3D), or 3-cell (element in 3D).

Return type:

None

See also

DMStagSetDOF

Source code at petsc4py/PETSc/DMStag.pyx:224

setGlobalSizes(sizes)#

Set global element counts in each dimension.

Logically collective.

Parameters:

sizes (tuple[int, ...]) – Global elementwise size in the one/two/three dimensions.

Return type:

None

Source code at petsc4py/PETSc/DMStag.pyx:246

setOwnershipRanges(ranges)#

Set elements per process in each dimension.

Logically collective.

Parameters:

ranges (tuple[Sequence[int], ...]) – Element counts for each process in one/two/three dimensions.

Return type:

None

Source code at petsc4py/PETSc/DMStag.pyx:286

setProcSizes(sizes)#

Set the number of processes in each dimension in the global process grid.

Logically collective.

Parameters:

sizes (tuple[int, ...]) – Number of processes in one/two/three dimensions.

Return type:

None

Source code at petsc4py/PETSc/DMStag.pyx:266

setStencilType(stenciltype)#

Set elementwise ghost/halo stencil type.

Logically collective.

Parameters:

stenciltype (StencilType | str) – The elementwise ghost stencil type.

Return type:

None

Source code at petsc4py/PETSc/DMStag.pyx:183

setStencilWidth(swidth)#

Set elementwise stencil width.

Logically collective.

The width value is not used when StencilType.NONE is specified.

Parameters:

swidth (int) – Stencil/halo/ghost width in elements.

Return type:

None

Source code at petsc4py/PETSc/DMStag.pyx:163

setUniformCoordinates(xmin=0, xmax=1, ymin=0, ymax=1, zmin=0, zmax=1)#

Set the coordinates to be a uniform grid..

Collective.

Local coordinates are populated, linearly extrapolated to ghost cells, including those outside the physical domain. This is also done in case of periodic boundaries, meaning that the same global point may have different coordinates in different local representations, which are equivalent assuming a periodicity implied by the arguments to this function, i.e., two points are equivalent if their difference is a multiple of xmax-xmin in the x dimension, ymax-ymin in the y dimension, and zmax-zmin in the z dimension.

Parameters:
  • xmin (float) – The minimum global coordinate value in the x dimension.

  • xmax (float) – The maximum global coordinate value in the x dimension.

  • ymin (float) – The minimum global coordinate value in the y dimension.

  • ymax (float) – The maximum global coordinate value in the y dimension.

  • zmin (float) – The minimum global coordinate value in the z dimension.

  • zmax (float) – The maximum global coordinate value in the z dimension.

Return type:

None

Source code at petsc4py/PETSc/DMStag.pyx:610

setUniformCoordinatesExplicit(xmin=0, xmax=1, ymin=0, ymax=1, zmin=0, zmax=1)#

Set coordinates to be a uniform grid, storing all values.

Collective.

Parameters:
  • xmin (float) – The minimum global coordinate value in the x dimension.

  • xmax (float) – The maximum global coordinate value in the x dimension.

  • ymin (float) – The minimum global coordinate value in the y dimension.

  • ymax (float) – The maximum global coordinate value in the y dimension.

  • zmin (float) – The minimum global coordinate value in the z dimension.

  • zmax (float) – The maximum global coordinate value in the z dimension.

Return type:

None

Source code at petsc4py/PETSc/DMStag.pyx:530

setUniformCoordinatesProduct(xmin=0, xmax=1, ymin=0, ymax=1, zmin=0, zmax=1)#

Create uniform coordinates, as a product of 1D arrays.

Collective.

The per-dimension 1-dimensional DMStag objects that comprise the product always have active 0-cells (vertices, element boundaries) and 1-cells (element centers).

Parameters:
  • xmin (float) – The minimum global coordinate value in the x dimension.

  • xmax (float) – The maximum global coordinate value in the x dimension.

  • ymin (float) – The minimum global coordinate value in the y dimension.

  • ymax (float) – The maximum global coordinate value in the y dimension.

  • zmin (float) – The minimum global coordinate value in the z dimension.

  • zmax (float) – The maximum global coordinate value in the z dimension.

Return type:

None

Source code at petsc4py/PETSc/DMStag.pyx:568

Attributes Documentation

boundary_types#

Boundary types in each dimension.

Source code at petsc4py/PETSc/DMStag.pyx:866

corners#

The lower left corner and size of local region in each dimension.

Source code at petsc4py/PETSc/DMStag.pyx:881

dim#

The dimension.

Source code at petsc4py/PETSc/DMStag.pyx:836

dofs#

The number of DOFs associated with each stratum of the grid.

Source code at petsc4py/PETSc/DMStag.pyx:841

entries_per_element#

The number of entries per element in the local representation.

Source code at petsc4py/PETSc/DMStag.pyx:846

ghost_corners#

The lower left corner and size of local region in each dimension.

Source code at petsc4py/PETSc/DMStag.pyx:886

global_sizes#

Global element counts in each dimension.

Source code at petsc4py/PETSc/DMStag.pyx:851

local_sizes#

Local elementwise sizes in each dimension.

Source code at petsc4py/PETSc/DMStag.pyx:856

proc_sizes#

The number of processes in each dimension in the global decomposition.

Source code at petsc4py/PETSc/DMStag.pyx:861

stencil_type#

Stencil type.

Source code at petsc4py/PETSc/DMStag.pyx:871

stencil_width#

Elementwise stencil width.

Source code at petsc4py/PETSc/DMStag.pyx:876