petsc4py.PETSc.DMDA#

class petsc4py.PETSc.DMDA#

Bases: DM

A DM object that is used to manage data for a structured grid.

Enumerations

ElementType

InterpolationType

StencilType

Methods Summary

create([dim, dof, sizes, proc_sizes, ...])

Create a DMDA object.

createNaturalVec()

Create a vector that will hold values in the natural numbering.

duplicate([dof, boundary_type, ...])

Duplicate a DMDA.

getAO()

Return the application ordering context for a distributed array.

getBoundaryType()

Return the type of ghost nodes at boundary in each dimension.

getCoordinateName(index)

Return the name of a coordinate dimension.

getCorners()

Return the lower left corner and the sizes of the owned local region.

getDim()

Return the topological dimension.

getDof()

Return the number of degrees of freedom per node.

getElementType()

Return the element type to be returned by getElements.

getElements([elem_type])

Return an array containing the indices of all the local elements.

getFieldName(field)

Return the name of an individual field component.

getGhostCorners()

Return the lower left corner and the size of the ghosted local region.

getGhostRanges()

Return the ranges of the local region in each dimension, including ghost nodes.

getInterpolationType()

Return the type of interpolation.

getOwnershipRanges()

Return the ranges of indices in each dimension owned by each process.

getProcSizes()

Return the number of processes in each dimension.

getRanges()

Return the ranges of the owned local region in each dimension.

getRefinementFactor()

Return the ratios that the DMDA grid is refined in each dimension.

getScatter()

Return the global-to-local, and local-to-local scatter contexts.

getSizes()

Return the number of grid points in each dimension.

getStencil()

Return the stencil type and width.

getStencilType()

Return the stencil type.

getStencilWidth()

Return the stencil width.

getVecArray(vec)

Get access to the vector.

globalToNatural(vg, vn[, addv])

Map values to the "natural" grid ordering.

naturalToGlobal(vn, vg[, addv])

Map values the to grid ordering.

setBoundaryType(boundary_type)

Set the type of ghost nodes on domain boundaries.

setCoordinateName(index, name)

Set the name of the coordinate dimension.

setDim(dim)

Set the topological dimension.

setDof(dof)

Set the number of degrees of freedom per vertex.

setElementType(elem_type)

Set the element type to be returned by getElements.

setFieldName(field, name)

Set the name of individual field components.

setInterpolationType(interp_type)

Set the type of interpolation.

setProcSizes(proc_sizes)

Set the number of processes in each dimension.

setRefinementFactor([refine_x, refine_y, ...])

Set the ratios for the DMDA grid refinement.

setSizes(sizes)

Set the number of grid points in each dimension.

setStencil(stencil_type, stencil_width)

Set the stencil type and width.

setStencilType(stencil_type)

Set the stencil type.

setStencilWidth(stencil_width)

Set the stencil width.

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

Set the DMDA coordinates to be a uniform grid.

Attributes Summary

boundary_type

Boundary types in each dimension.

corners

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

dim

The grid dimension.

dof

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

ghost_corners

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

ghost_ranges

Ranges of local region, including ghost nodes.

proc_sizes

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

ranges

Ranges of the local region in each dimension.

sizes

The global dimension.

stencil

Stencil type and width.

stencil_type

Stencil type.

stencil_width

Elementwise stencil width.

Methods Documentation

create(dim=None, dof=None, sizes=None, proc_sizes=None, boundary_type=None, stencil_type=None, stencil_width=None, setup=True, ownership_ranges=None, comm=None)#

Create a DMDA object.

Collective.

This routine performs the following steps of the C API: - petsc.DMDACreate - petsc.DMSetDimension - petsc.DMDASetDof - petsc.DMDASetSizes - petsc.DMDASetNumProcs - petsc.DMDASetOwnershipRanges - petsc.DMDASetBoundaryType - petsc.DMDASetStencilType - petsc.DMDASetStencilWidth - petsc.DMSetUp (optionally)

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

  • dof (int | None) – The number of degrees of freedom.

  • sizes (DimsSpec | None) – The number of elements in each dimension.

  • proc_sizes (DimsSpec | None) – The number of processes in x, y, z dimensions.

  • boundary_type (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.

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

  • 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.

Return type:

Self

Source code at petsc4py/PETSc/DMDA.pyx:26

createNaturalVec()#

Create a vector that will hold values in the natural numbering.

Collective.

The number of local entries in the vector on each process is the same as in a vector created with DM.createGlobalVec.

Source code at petsc4py/PETSc/DMDA.pyx:827

Return type:

Vec

duplicate(dof=None, boundary_type=None, stencil_type=None, stencil_width=None)#

Duplicate a DMDA.

Collective.

This routine retrieves the information from the DMDA and recreates it. Parameters dof, boundary_type, stencil_type, stencil_width will be overwritten, if provided.

Parameters:
Return type:

DMDA

Source code at petsc4py/PETSc/DMDA.pyx:143

getAO()#

Return the application ordering context for a distributed array.

Collective.

The returned AO maps to the natural grid ordering that would be used for the DMDA if only 1 processor were employed (ordering most rapidly in the x-dimension, then y, then z). Multiple degrees of freedom are numbered for each node (rather than 1 component for the whole grid, then the next component, etc.).

See also

DMDAGetAO

Source code at petsc4py/PETSc/DMDA.pyx:906

Return type:

AO

getBoundaryType()#

Return the type of ghost nodes at boundary in each dimension.

Not collective.

See also

setBoundaryType

Source code at petsc4py/PETSc/DMDA.pyx:404

Return type:

tuple[BoundaryType, …]

getCoordinateName(index)#

Return the name of a coordinate dimension.

Not collective.

Parameters:

index (int) – The coordinate number for the DMDA (0, 1, …, dim-1).

Return type:

str

Source code at petsc4py/PETSc/DMDA.pyx:805

getCorners()#

Return the lower left corner and the sizes of the owned local region.

Not collective.

Returns the global (x,y,z) indices of the lower left corner (first tuple) and size of the local region (second tuple).

Excluding ghost points.

The corner information is independent of the number of degrees of freedom per node. Thus the returned values can be thought of as coordinates on a logical grid, where each grid point has (potentially) several degrees of freedom.

Source code at petsc4py/PETSc/DMDA.pyx:620

Return type:

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

getDim()#

Return the topological dimension.

Not collective.

Source code at petsc4py/PETSc/DMDA.pyx:230

Return type:

int

getDof()#

Return the number of degrees of freedom per node.

Not collective.

See also

setDof, DMDAGetInfo

Source code at petsc4py/PETSc/DMDA.pyx:260

Return type:

int

getElementType()#

Return the element type to be returned by getElements.

Not collective.

Source code at petsc4py/PETSc/DMDA.pyx:1048

Return type:

ElementType

getElements(elem_type=None)#

Return an array containing the indices of all the local elements.

Not collective.

The elements are in local coordinates.

Each process uniquely owns a subset of the elements. That is, no element is owned by two or more processes.

Parameters:

elem_type (ElementType | None) – The element type.

Return type:

ArrayInt

See also

DMDAGetElements

Source code at petsc4py/PETSc/DMDA.pyx:1062

getFieldName(field)#

Return the name of an individual field component.

Not collective.

Parameters:

field (int) – The field number for the DMDA (0, 1, …, dof-1), where dof indicates the number of degrees of freedom per node within the DMDA.

Return type:

str

Source code at petsc4py/PETSc/DMDA.pyx:697

getGhostCorners()#

Return the lower left corner and the size of the ghosted local region.

Not collective.

Returns the global (x,y,z) indices of the lower left corner (first tuple) and size of the local region (second tuple).

Source code at petsc4py/PETSc/DMDA.pyx:649

Return type:

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

getGhostRanges()#

Return the ranges of the local region in each dimension, including ghost nodes.

Not collective.

Source code at petsc4py/PETSc/DMDA.pyx:575

Return type:

tuple[tuple[int, int], …]

getInterpolationType()#

Return the type of interpolation.

Not collective.

Source code at petsc4py/PETSc/DMDA.pyx:1018

Return type:

InterpolationType

getOwnershipRanges()#

Return the ranges of indices in each dimension owned by each process.

Not collective.

These numbers are not multiplied by the number of DOFs per node.

Source code at petsc4py/PETSc/DMDA.pyx:595

Return type:

tuple[ArrayInt, …]

getProcSizes()#

Return the number of processes in each dimension.

Not collective.

Source code at petsc4py/PETSc/DMDA.pyx:357

Return type:

tuple[int, …]

getRanges()#

Return the ranges of the owned local region in each dimension.

Not collective.

Excluding ghost nodes.

Source code at petsc4py/PETSc/DMDA.pyx:553

Return type:

tuple[tuple[int, int], …]

getRefinementFactor()#

Return the ratios that the DMDA grid is refined in each dimension.

Not collective.

Source code at petsc4py/PETSc/DMDA.pyx:979

Return type:

tuple[int, …]

getScatter()#

Return the global-to-local, and local-to-local scatter contexts.

Collective.

See also

DMDAGetScatter

Source code at petsc4py/PETSc/DMDA.pyx:927

Return type:

tuple[Scatter, Scatter]

getSizes()#

Return the number of grid points in each dimension.

Not collective.

See also

setSizes, DMDAGetInfo

Source code at petsc4py/PETSc/DMDA.pyx:307

Return type:

tuple[int, …]

getStencil()#

Return the stencil type and width.

Not collective.

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

Return type:

tuple[StencilType, int]

getStencilType()#

Return the stencil type.

Not collective.

Source code at petsc4py/PETSc/DMDA.pyx:445

Return type:

StencilType

getStencilWidth()#

Return the stencil width.

Not collective.

See also

setStencilWidth

Source code at petsc4py/PETSc/DMDA.pyx:483

Return type:

int

getVecArray(vec)#

Get access to the vector.

Not collective.

Use via The with statement context manager (PEP 343).

Parameters:

vec (Vec) – The vector to which access is being requested.

Return type:

Any

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

globalToNatural(vg, vn, addv=None)#

Map values to the “natural” grid ordering.

Neighborwise collective.

You must call createNaturalVec before using this routine.

Parameters:
  • vg (Vec) – The global vector in a grid ordering.

  • vn (Vec) – The global vector in a “natural” ordering.

  • addv (InsertMode | None) – The insertion mode.

Return type:

None

Source code at petsc4py/PETSc/DMDA.pyx:844

naturalToGlobal(vn, vg, addv=None)#

Map values the to grid ordering.

Neighborwise collective.

Parameters:
  • vn (Vec) – The global vector in a natural ordering.

  • vg (Vec) – the global vector in a grid ordering.

  • addv (InsertMode | None) – The insertion mode.

Return type:

None

Source code at petsc4py/PETSc/DMDA.pyx:875

setBoundaryType(boundary_type)#

Set the type of ghost nodes on domain boundaries.

Not collective.

Parameters:

boundary_type (tuple[BoundaryType | int | str | bool, ...]) – The boundary type in (x), (x, y), or (x, y, z) dimensions.

Return type:

None

Source code at petsc4py/PETSc/DMDA.pyx:380

setCoordinateName(index, name)#

Set the name of the coordinate dimension.

Logically collective.

Parameters:
  • index (int) – The coordinate number for the DMDA (0, 1, …, dim-1).

  • name (str) – The name of the coordinate.

Return type:

None

Source code at petsc4py/PETSc/DMDA.pyx:783

setDim(dim)#

Set the topological dimension.

Collective.

Parameters:

dim (int) – Topological dimension.

Return type:

None

Source code at petsc4py/PETSc/DMDA.pyx:213

setDof(dof)#

Set the number of degrees of freedom per vertex.

Not collective.

Parameters:

dof (int) – The number of degrees of freedom.

Return type:

None

See also

getDof, DMDASetDof

Source code at petsc4py/PETSc/DMDA.pyx:242

setElementType(elem_type)#

Set the element type to be returned by getElements.

Not collective.

Source code at petsc4py/PETSc/DMDA.pyx:1034

Parameters:

elem_type (ElementType | str) –

Return type:

None

setFieldName(field, name)#

Set the name of individual field components.

Logically collective.

Parameters:
  • field (int) – The field number for the DMDA (0, 1, …, dof-1), where dof indicates the number of degrees of freedom per node within the DMDA.

  • name (str) – The name of the field (component).

Return type:

None

Source code at petsc4py/PETSc/DMDA.pyx:673

setInterpolationType(interp_type)#

Set the type of interpolation.

Logically collective.

You should call this on the coarser of the two DMDAs you pass to DM.createInterpolation.

Parameters:

interp_type (InterpolationType) – The interpolation type.

Return type:

None

Source code at petsc4py/PETSc/DMDA.pyx:997

setProcSizes(proc_sizes)#

Set the number of processes in each dimension.

Logically collective.

Parameters:

proc_sizes (DimsSpec) – The number of processes in (x,), (x, y), or (x, y, z) dimensions.

Return type:

None

Source code at petsc4py/PETSc/DMDA.pyx:330

setRefinementFactor(refine_x=2, refine_y=2, refine_z=2)#

Set the ratios for the DMDA grid refinement.

Logically collective.

Parameters:
  • refine_x (int) – Ratio of fine grid to coarse in x dimension.

  • refine_y (int) – Ratio of fine grid to coarse in y dimension.

  • refine_z (int) – Ratio of fine grid to coarse in z dimension.

Return type:

None

Source code at petsc4py/PETSc/DMDA.pyx:946

setSizes(sizes)#

Set the number of grid points in each dimension.

Logically collective.

Parameters:

sizes (DimsSpec) – The global (x,), (x, y), or (x, y, z) size.

Return type:

None

Source code at petsc4py/PETSc/DMDA.pyx:280

setStencil(stencil_type, stencil_width)#

Set the stencil type and width.

Not collective.

Parameters:
  • stencil_type (StencilType) – The stencil type.

  • stencil_width (int) – The stencil width.

Return type:

None

Source code at petsc4py/PETSc/DMDA.pyx:503

setStencilType(stencil_type)#

Set the stencil type.

Logically collective.

Parameters:
  • stype – The stencil type.

  • stencil_type (StencilType) –

Return type:

None

Source code at petsc4py/PETSc/DMDA.pyx:427

setStencilWidth(stencil_width)#

Set the stencil width.

Logically collective.

Parameters:

stencil_width (int) – The stencil width.

Return type:

None

Source code at petsc4py/PETSc/DMDA.pyx:465

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

Set the DMDA coordinates to be a uniform grid.

Collective.

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

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

  • ymin (float) – The minimum in the y dimension (value ignored for 1 dimensional problems).

  • ymax (float) – The maximum in the y dimension (value ignored for 1 dimensional problems).

  • zmin (float) – The minimum in the z dimension (value ignored for 1 or 2 dimensional problems).

  • zmax (float) – The maximum in the z dimension (value ignored for 1 or 2 dimensional problems).

Return type:

None

Source code at petsc4py/PETSc/DMDA.pyx:738

Attributes Documentation

boundary_type#

Boundary types in each dimension.

Source code at petsc4py/PETSc/DMDA.pyx:1121

corners#

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

Source code at petsc4py/PETSc/DMDA.pyx:1151

dim#

The grid dimension.

Source code at petsc4py/PETSc/DMDA.pyx:1101

dof#

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

Source code at petsc4py/PETSc/DMDA.pyx:1106

ghost_corners#

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

Source code at petsc4py/PETSc/DMDA.pyx:1156

ghost_ranges#

Ranges of local region, including ghost nodes.

Source code at petsc4py/PETSc/DMDA.pyx:1146

proc_sizes#

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

Source code at petsc4py/PETSc/DMDA.pyx:1116

ranges#

Ranges of the local region in each dimension.

Source code at petsc4py/PETSc/DMDA.pyx:1141

sizes#

The global dimension.

Source code at petsc4py/PETSc/DMDA.pyx:1111

stencil#

Stencil type and width.

Source code at petsc4py/PETSc/DMDA.pyx:1126

stencil_type#

Stencil type.

Source code at petsc4py/PETSc/DMDA.pyx:1131

stencil_width#

Elementwise stencil width.

Source code at petsc4py/PETSc/DMDA.pyx:1136