petsc4py.PETSc.DM#

class petsc4py.PETSc.DM#

Bases: Object

An object describing a computational grid or mesh.

Enumerations

BoundaryType

DM Boundary types.

PolytopeType

The DM cell types.

ReorderDefaultFlag

The DM reordering default flags.

Type

DM types.

Methods Summary

adaptLabel(label)

Adapt a DM based on a DMLabel.

adaptMetric(metric[, bdLabel, rgLabel])

Return a mesh adapted to the specified metric field.

addCoarsenHook(coarsenhook, restricthook[, ...])

Add a callback to be executed when restricting to a coarser grid.

addField(field[, label])

Add a field to a DM object.

appendOptionsPrefix(prefix)

Append to the prefix used for searching for options in the database.

clearDS()

Remove all discrete systems from the DM.

clearFields()

Remove all fields from the DM.

clearLabelStratum(name, value)

Remove all points from a stratum.

clearLabelValue(name, point, value)

Remove a point from a DMLabel with given value.

clone()

Return the cloned DM .

coarsen([comm])

Return a coarsened DM object.

coarsenHierarchy(nlevels)

Coarsen this DM and return the coarsened DM hierarchy.

convert(dm_type)

Return a DM converted to another DM.

copyDS(dm)

Copy the discrete systems for this DM into another DM.

copyDisc(dm)

Copy fields and discrete systems of a DM into another DM.

copyFields(dm)

Copy the discretizations of this DM into another DM.

create([comm])

Return an empty DM.

createDS()

Create discrete systems.

createFieldDecomposition()

Return a list of IS objects.

createGlobalVec()

Return a global vector.

createInjection(dm)

Return the injection matrix into a finer DM.

createInterpolation(dm)

Return the interpolation matrix to a finer DM.

createLabel(name)

Create a label of the given name if it does not already exist.

createLocalVec()

Return a local vector.

createMassMatrix(dmf)

Return the mass matrix between this DM and the given DM.

createMat()

Return an empty matrix.

createRestriction(dm)

Return the restriction matrix between this DM and the given DM.

createSectionSF(localsec, globalsec)

Create the SF encoding the parallel DOF overlap for the DM.

createSubDM(fields)

Return IS and DM encapsulating a subproblem.

destroy()

Destroy the object.

getAppCtx()

Return the application context.

getAuxiliaryVec([label, value, part])

Return an auxiliary vector for region.

getBasicAdjacency()

Return the flags for determining variable influence.

getBlockSize()

Return the inherent block size associated with a DM.

getBoundingBox()

Return the dimension of embedding space for coordinates values.

getCellCoordinateDM()

Return the cell coordinate DM.

getCellCoordinateSection()

Return the cell coordinate layout over the DM.

getCellCoordinates()

Return a global vector with the cellwise coordinates.

getCellCoordinatesLocal()

Return a local vector with the cellwise coordinates.

getCoarseDM()

Return the coarse DM.

getCoarsenLevel()

Return the number of coarsenings.

getCoordinateDM()

Return the coordinate DM.

getCoordinateDim()

Return the dimension of embedding space for coordinates values.

getCoordinateSection()

Return coordinate values layout over the mesh.

getCoordinates()

Return a global vector with the coordinates associated.

getCoordinatesLocal()

Return a local vector with the coordinates associated.

getCoordinatesLocalized()

Check if the coordinates have been localized for cells.

getDS()

Return default DS.

getDimension()

Return the topological dimension of the DM.

getField(index)

Return the discretization object for a given DM field.

getFieldAdjacency(field)

Return the flags for determining variable influence.

getGlobalSection()

Return the Section encoding the global data layout for the DM.

getGlobalVec()

Return a global vector.

getLGMap()

Return local mapping to global mapping.

getLabel(name)

Return the label of a given name.

getLabelIdIS(name)

Return an IS of all values that the DMLabel takes.

getLabelName(index)

Return the name of nth label.

getLabelOutput(name)

Return the output flag for a given label.

getLabelSize(name)

Return the number of values that the DMLabel takes.

getLabelValue(name, point)

Return the value in DMLabel for the given point.

getLocalBoundingBox()

Return the bounding box for the piece of the DM.

getLocalSection()

Return the Section encoding the local data layout for the DM.

getLocalVec()

Return a local vector.

getNumFields()

Return the number of fields in the DM.

getNumLabels()

Return the number of labels defined by on the DM.

getOptionsPrefix()

Return the prefix used for searching for options in the database.

getPointSF()

Return the SF encoding the parallel DOF overlap for the DM.

getRefineLevel()

Return the refinement level.

getSectionSF()

Return the Section encoding the parallel DOF overlap.

getStratumIS(name, value)

Return the points in a label stratum.

getStratumSize(name, value)

Return the number of points in a label stratum.

getType()

Return the DM type name.

globalToLocal(vg, vl[, addv])

Update local vectors from global vector.

hasLabel(name)

Determine whether the DM has a label.

load(viewer)

Return a DM stored in binary.

localToGlobal(vl, vg[, addv])

Update global vectors from local vector.

localToLocal(vl, vlg[, addv])

Map the values from a local vector to another local vector.

localizeCoordinates()

Create local coordinates for cells having periodic faces.

refine([comm])

Return a refined DM object.

refineHierarchy(nlevels)

Refine this DM and return the refined DM hierarchy.

removeLabel(name)

Remove and destroy the label by name.

restoreGlobalVec(vg)

Restore a global vector.

restoreLocalVec(vl)

Restore a local vector.

setAppCtx(appctx)

Set the application context.

setAuxiliaryVec(aux, label[, value, part])

Set an auxiliary vector for a specific region.

setBasicAdjacency(useCone, useClosure)

Set the flags for determining variable influence.

setCellCoordinateDM(dm)

Set the cell coordinate DM.

setCellCoordinateSection(dim, sec)

Set the cell coordinate layout over the DM.

setCellCoordinates(c)

Set a global vector with the cellwise coordinates.

setCellCoordinatesLocal(c)

Set a local vector with the cellwise coordinates.

setCoarseDM(dm)

Set the coarse DM.

setCoordinateDim(dim)

Set the dimension of embedding space for coordinates values.

setCoordinateDisc(disc, project)

Project coordinates to a different space.

setCoordinates(c)

Set a global vector that holds the coordinates.

setCoordinatesLocal(c)

Set a local vector with the ghost point holding the coordinates.

setDimension(dim)

Set the topological dimension of the DM.

setField(index, field[, label])

Set the discretization object for a given DM field.

setFieldAdjacency(field, useCone, useClosure)

Set the flags for determining variable influence.

setFromOptions()

Configure the object from the options database.

setGlobalSection(sec)

Set the Section encoding the global data layout for the DM.

setKSPComputeOperators(operators[, args, kargs])

Matrix associated with the linear system.

setLabelOutput(name, output)

Set if a given label should be saved to a view.

setLabelValue(name, point, value)

Set a point to a DMLabel with a given value.

setLocalSection(sec)

Set the Section encoding the local data layout for the DM.

setMatType(mat_type)

Set matrix type to be used by DM.createMat.

setNumFields(numFields)

Set the number of fields in the DM.

setOptionsPrefix(prefix)

Set the prefix used for searching for options in the database.

setPointSF(sf)

Set the SF encoding the parallel DOF overlap for the DM.

setRefineLevel(level)

Set the number of refinements.

setSNESFunction(function[, args, kargs])

Set SNES residual evaluation function.

setSNESJacobian(jacobian[, args, kargs])

Set the SNES Jacobian evaluation function.

setSectionSF(sf)

Set the Section encoding the parallel DOF overlap for the DM.

setType(dm_type)

Build a DM.

setUp()

Return the data structure.

setVecType(vec_type)

Set the type of vector.

view([viewer])

View the DM.

Attributes Summary

appctx

Application context.

ds

Discrete space.

Methods Documentation

adaptLabel(label)#

Adapt a DM based on a DMLabel.

Collective.

Parameters:

label (str) – The name of the DMLabel.

Return type:

DM

See also

DMAdaptLabel

Source code at petsc4py/PETSc/DM.pyx:1571

adaptMetric(metric, bdLabel=None, rgLabel=None)#

Return a mesh adapted to the specified metric field.

Collective.

Parameters:
  • metric (Vec) – The metric to which the mesh is adapted, defined vertex-wise.

  • bdLabel (str | None) – Label for boundary tags.

  • rgLabel (str | None) – Label for cell tag.

Return type:

DM

See also

DMAdaptMetric

Source code at petsc4py/PETSc/DM.pyx:1594

addCoarsenHook(coarsenhook, restricthook, args=None, kargs=None)#

Add a callback to be executed when restricting to a coarser grid.

Logically collective.

Parameters:
Return type:

None

See also

DMCoarsenHookAdd

Source code at petsc4py/PETSc/DM.pyx:2259

addField(field, label=None)#

Add a field to a DM object.

Logically collective.

Parameters:
  • field (Object) – The discretization object.

  • label (str | None) – The name of the label indicating the support of the field, or None for the entire mesh.

Return type:

None

See also

DMAddField

Source code at petsc4py/PETSc/DM.pyx:611

appendOptionsPrefix(prefix)#

Append to the prefix used for searching for options in the database.

Logically collective.

Source code at petsc4py/PETSc/DM.pyx:298

Parameters:

prefix (str | None)

Return type:

None

clearDS()#

Remove all discrete systems from the DM.

Logically collective.

See also

DMClearDS

Source code at petsc4py/PETSc/DM.pyx:675

Return type:

None

clearFields()#

Remove all fields from the DM.

Logically collective.

See also

DMClearFields

Source code at petsc4py/PETSc/DM.pyx:634

Return type:

None

clearLabelStratum(name, value)#

Remove all points from a stratum.

Not collective.

Parameters:
  • name (str) – The label name.

  • value (int) – The stratum value.

Return type:

None

Source code at petsc4py/PETSc/DM.pyx:2050

clearLabelValue(name, point, value)#

Remove a point from a DMLabel with given value.

Not collective.

Parameters:
  • name (str) – The label name.

  • point (int) – The mesh point.

  • value (int) – The label value for the point.

Return type:

None

Source code at petsc4py/PETSc/DM.pyx:1936

clone()#

Return the cloned DM .

Collective.

See also

DMClone

Source code at petsc4py/PETSc/DM.pyx:155

Return type:

DM

coarsen(comm=None)#

Return a coarsened DM object.

Collective.

Parameters:

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

Return type:

DM

See also

DMCoarsen

Source code at petsc4py/PETSc/DM.pyx:1445

coarsenHierarchy(nlevels)#

Coarsen this DM and return the coarsened DM hierarchy.

Collective.

Parameters:

nlevels (int) – The number of levels of coarsening.

Return type:

list

Source code at petsc4py/PETSc/DM.pyx:1496

convert(dm_type)#

Return a DM converted to another DM.

Collective.

Parameters:

dm_type (Type | str) – The new DM.Type, use “same” for the same type.

Return type:

DM

See also

DM.Type, DMConvert

Source code at petsc4py/PETSc/DM.pyx:1398

copyDS(dm)#

Copy the discrete systems for this DM into another DM.

Collective.

Parameters:

dm (DM) – The DM that the discrete fields are copied into.

Return type:

None

See also

DMCopyDS

Source code at petsc4py/PETSc/DM.pyx:702

copyDisc(dm)#

Copy fields and discrete systems of a DM into another DM.

Collective.

Parameters:

dm (DM) – The DM that the fields and discrete systems are copied into.

Return type:

None

See also

DMCopyDisc

Source code at petsc4py/PETSc/DM.pyx:719

copyFields(dm)#

Copy the discretizations of this DM into another DM.

Collective.

Parameters:

dm (DM) – The DM that the fields are copied into.

Return type:

None

See also

DMCopyFields

Source code at petsc4py/PETSc/DM.pyx:646

create(comm=None)#

Return an empty DM.

Collective.

Parameters:

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

Return type:

Self

See also

DMCreate

Source code at petsc4py/PETSc/DM.pyx:134

createDS()#

Create discrete systems.

Collective.

See also

DMCreateDS

Source code at petsc4py/PETSc/DM.pyx:663

Return type:

None

createFieldDecomposition()#

Return a list of IS objects.

Not collective.

Notes

The user is responsible for freeing all requested arrays.

Source code at petsc4py/PETSc/DM.pyx:2148

Return type:

tuple[list, list, list]

createGlobalVec()#

Return a global vector.

Collective.

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

Return type:

Vec

createInjection(dm)#

Return the injection matrix into a finer DM.

Collective.

Parameters:

dm (DM) – The second, finer DM object.

Return type:

Mat

Source code at petsc4py/PETSc/DM.pyx:1360

createInterpolation(dm)#

Return the interpolation matrix to a finer DM.

Collective.

Parameters:

dm (DM) – The second, finer DM.

Return type:

tuple[Mat, Vec]

Source code at petsc4py/PETSc/DM.pyx:1339

createLabel(name)#

Create a label of the given name if it does not already exist.

Not collective.

Parameters:

name (str) – The label name.

Return type:

None

See also

DMCreateLabel

Source code at petsc4py/PETSc/DM.pyx:1848

createLocalVec()#

Return a local vector.

Not collective.

Source code at petsc4py/PETSc/DM.pyx:780

Return type:

Vec

createMassMatrix(dmf)#

Return the mass matrix between this DM and the given DM.

Collective.

Parameters:

dmf (DM) – The second DM.

Return type:

Mat

Source code at petsc4py/PETSc/DM.pyx:1320

createMat()#

Return an empty matrix.

Collective.

See also

DMCreateMatrix

Source code at petsc4py/PETSc/DM.pyx:1306

Return type:

Mat

createRestriction(dm)#

Return the restriction matrix between this DM and the given DM.

Collective.

Parameters:

dm (DM) – The second, finer DM object.

Return type:

Mat

Source code at petsc4py/PETSc/DM.pyx:1379

createSectionSF(localsec, globalsec)#

Create the SF encoding the parallel DOF overlap for the DM.

Collective.

Parameters:
  • localsec (Section) – Describe the local data layout.

  • globalsec (Section) – Describe the global data layout.

Return type:

None

Notes

Encoding based on the Section describing the data layout.

Source code at petsc4py/PETSc/DM.pyx:1712

createSubDM(fields)#

Return IS and DM encapsulating a subproblem.

Not collective.

Returns:

  • iset (IS) – The global indices for all the degrees of freedom.

  • subdm (DM) – The DM for the subproblem.

Parameters:

fields (Sequence[int])

Return type:

tuple[IS, DM]

See also

DMCreateSubDM

Source code at petsc4py/PETSc/DM.pyx:446

destroy()#

Destroy the object.

Collective.

See also

DMDestroy

Source code at petsc4py/PETSc/DM.pyx:121

Return type:

Self

getAppCtx()#

Return the application context.

Source code at petsc4py/PETSc/DM.pyx:343

Return type:

Any

getAuxiliaryVec(label=None, value=0, part=0)#

Return an auxiliary vector for region.

Not collective.

Parameters:
  • label (str | None) – The name of the DMLabel.

  • value (int | None) – Indicate the region.

  • part (int | None) – The equation part, or 0 is unused.

Return type:

Vec

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

getBasicAdjacency()#

Return the flags for determining variable influence.

Not collective.

Returns:

  • useCone (bool) – Whether adjacency uses cone information.

  • useClosure (bool) – Whether adjacency is computed using full closure information.

Return type:

tuple[bool, bool]

Source code at petsc4py/PETSc/DM.pyx:370

getBlockSize()#

Return the inherent block size associated with a DM.

Not collective.

See also

DMGetBlockSize

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

Return type:

int

getBoundingBox()#

Return the dimension of embedding space for coordinates values.

Not collective.

See also

DMGetBoundingBox

Source code at petsc4py/PETSc/DM.pyx:1232

Return type:

tuple[tuple[float, float], …]

getCellCoordinateDM()#

Return the cell coordinate DM.

Collective.

Source code at petsc4py/PETSc/DM.pyx:1085

Return type:

DM

getCellCoordinateSection()#

Return the cell coordinate layout over the DM.

Collective.

Source code at petsc4py/PETSc/DM.pyx:1120

Return type:

Section

getCellCoordinates()#

Return a global vector with the cellwise coordinates.

Collective.

Source code at petsc4py/PETSc/DM.pyx:1152

Return type:

Vec

getCellCoordinatesLocal()#

Return a local vector with the cellwise coordinates.

Collective.

Source code at petsc4py/PETSc/DM.pyx:1184

Return type:

Vec

getCoarseDM()#

Return the coarse DM.

Collective.

See also

DMGetCoarseDM

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

Return type:

DM

getCoarsenLevel()#

Return the number of coarsenings.

Not collective.

Source code at petsc4py/PETSc/DM.pyx:1555

Return type:

int

getCoordinateDM()#

Return the coordinate DM.

Collective.

Source code at petsc4py/PETSc/DM.pyx:974

Return type:

DM

getCoordinateDim()#

Return the dimension of embedding space for coordinates values.

Not collective.

Source code at petsc4py/PETSc/DM.pyx:238

Return type:

int

getCoordinateSection()#

Return coordinate values layout over the mesh.

Collective.

Source code at petsc4py/PETSc/DM.pyx:989

Return type:

Section

getCoordinates()#

Return a global vector with the coordinates associated.

Collective.

See also

DMGetCoordinates

Source code at petsc4py/PETSc/DM.pyx:1021

Return type:

Vec

getCoordinatesLocal()#

Return a local vector with the coordinates associated.

Collective the first time it is called.

Source code at petsc4py/PETSc/DM.pyx:1053

Return type:

Vec

getCoordinatesLocalized()#

Check if the coordinates have been localized for cells.

Not collective.

Source code at petsc4py/PETSc/DM.pyx:1218

Return type:

bool

getDS()#

Return default DS.

Not collective.

See also

DMGetDS

Source code at petsc4py/PETSc/DM.pyx:687

Return type:

DS

getDimension()#

Return the topological dimension of the DM.

Not collective.

See also

DMGetDimension

Source code at petsc4py/PETSc/DM.pyx:206

Return type:

int

getField(index)#

Return the discretization object for a given DM field.

Not collective.

Parameters:

index (int) – The field number.

Return type:

tuple[Object, None]

See also

DMGetField

Source code at petsc4py/PETSc/DM.pyx:586

getFieldAdjacency(field)#

Return the flags for determining variable influence.

Not collective.

Parameters:

field (int) – The field number.

Returns:

  • useCone (bool) – Whether adjacency uses cone information.

  • useClosure (bool) – Whether adjacency is computed using full closure information.

Return type:

tuple[bool, bool]

See also

DMGetAdjacency

Source code at petsc4py/PETSc/DM.pyx:416

getGlobalSection()#

Return the Section encoding the global data layout for the DM.

Collective the first time it is called.

Source code at petsc4py/PETSc/DM.pyx:1688

Return type:

Section

getGlobalVec()#

Return a global vector.

Collective.

Source code at petsc4py/PETSc/DM.pyx:794

Return type:

Vec

getLGMap()#

Return local mapping to global mapping.

Collective.

Source code at petsc4py/PETSc/DM.pyx:929

Return type:

LGMap

getLabel(name)#

Return the label of a given name.

Not collective.

See also

DMGetLabel

Source code at petsc4py/PETSc/DM.pyx:1630

Parameters:

name (str)

Return type:

DMLabel

getLabelIdIS(name)#

Return an IS of all values that the DMLabel takes.

Not collective.

Parameters:

name (str) – The label name.

Return type:

IS

Source code at petsc4py/PETSc/DM.pyx:1981

getLabelName(index)#

Return the name of nth label.

Not collective.

Parameters:

index (int) – The label number.

Return type:

str

See also

DMGetLabelName

Source code at petsc4py/PETSc/DM.pyx:1807

getLabelOutput(name)#

Return the output flag for a given label.

Not collective.

Parameters:

name (str) – The label name.

Return type:

bool

See also

DMGetLabelOutput

Source code at petsc4py/PETSc/DM.pyx:2094

getLabelSize(name)#

Return the number of values that the DMLabel takes.

Not collective.

Parameters:

name (str) – The label name.

Return type:

int

Source code at petsc4py/PETSc/DM.pyx:1960

getLabelValue(name, point)#

Return the value in DMLabel for the given point.

Not collective.

Parameters:
  • name (str) – The label name.

  • point (int) – The mesh point

Return type:

int

See also

DMGetLabelValue

Source code at petsc4py/PETSc/DM.pyx:1889

getLocalBoundingBox()#

Return the bounding box for the piece of the DM.

Not collective.

Source code at petsc4py/PETSc/DM.pyx:1249

Return type:

tuple[tuple[float, float], …]

getLocalSection()#

Return the Section encoding the local data layout for the DM.

Not collective.

Source code at petsc4py/PETSc/DM.pyx:1661

Return type:

Section

getLocalVec()#

Return a local vector.

Not collective.

See also

DMGetLocalVector

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

Return type:

Vec

getNumFields()#

Return the number of fields in the DM.

Not collective.

See also

DMGetNumFields

Source code at petsc4py/PETSc/DM.pyx:546

Return type:

int

getNumLabels()#

Return the number of labels defined by on the DM.

Not collective.

See also

DMGetNumLabels

Source code at petsc4py/PETSc/DM.pyx:1793

Return type:

int

getOptionsPrefix()#

Return the prefix used for searching for options in the database.

Not collective.

Source code at petsc4py/PETSc/DM.pyx:284

Return type:

str

getPointSF()#

Return the SF encoding the parallel DOF overlap for the DM.

Not collective.

See also

DMGetPointSF

Source code at petsc4py/PETSc/DM.pyx:1766

Return type:

SF

getRefineLevel()#

Return the refinement level.

Not collective.

See also

DMGetRefineLevel

Source code at petsc4py/PETSc/DM.pyx:1523

Return type:

int

getSectionSF()#

Return the Section encoding the parallel DOF overlap.

Collective the first time it is called.

See also

DMGetSectionSF

Source code at petsc4py/PETSc/DM.pyx:1735

Return type:

SF

getStratumIS(name, value)#

Return the points in a label stratum.

Not collective.

Parameters:
  • name (str) – The label name.

  • value (int) – The stratum value.

Return type:

IS

See also

DMGetStratumIS

Source code at petsc4py/PETSc/DM.pyx:2026

getStratumSize(name, value)#

Return the number of points in a label stratum.

Not collective.

Parameters:
  • name (str) – The label name.

  • value (int) – The stratum value.

Return type:

int

See also

DMGetStratumSize

Source code at petsc4py/PETSc/DM.pyx:2002

getType()#

Return the DM type name.

Not collective.

See also

DMGetType

Source code at petsc4py/PETSc/DM.pyx:192

Return type:

str

globalToLocal(vg, vl, addv=None)#

Update local vectors from global vector.

Neighborwise collective.

Parameters:
Return type:

None

Source code at petsc4py/PETSc/DM.pyx:860

hasLabel(name)#

Determine whether the DM has a label.

Not collective.

Parameters:

name (str) – The label name.

Return type:

bool

See also

DMHasLabel

Source code at petsc4py/PETSc/DM.pyx:1827

load(viewer)#

Return a DM stored in binary.

Collective.

Parameters:

viewer (Viewer) – Viewer used to store the DM, like Viewer.Type.BINARY or Viewer.Type.HDF5.

Return type:

Self

Notes

When using Viewer.Type.HDF5 format, one can save multiple DMPlex meshes in a single HDF5 files. This in turn requires one to name the DMPlex object with Object.setName before saving it with DM.view and before loading it with DM.load for identification of the mesh object.

Source code at petsc4py/PETSc/DM.pyx:95

localToGlobal(vl, vg, addv=None)#

Update global vectors from local vector.

Neighborwise collective.

Parameters:
Return type:

None

Source code at petsc4py/PETSc/DM.pyx:883

localToLocal(vl, vlg, addv=None)#

Map the values from a local vector to another local vector.

Neighborwise collective.

Parameters:
Return type:

None

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

localizeCoordinates()#

Create local coordinates for cells having periodic faces.

Collective.

Notes

Used if the mesh is periodic.

Source code at petsc4py/PETSc/DM.pyx:1266

Return type:

None

refine(comm=None)#

Return a refined DM object.

Collective.

Parameters:

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

Return type:

DM

See also

DMRefine

Source code at petsc4py/PETSc/DM.pyx:1421

refineHierarchy(nlevels)#

Refine this DM and return the refined DM hierarchy.

Collective.

Parameters:

nlevels (int) – The number of levels of refinement.

Return type:

list

Source code at petsc4py/PETSc/DM.pyx:1469

removeLabel(name)#

Remove and destroy the label by name.

Not collective.

Parameters:

name (str) – The label name.

Return type:

None

See also

DMRemoveLabel

Source code at petsc4py/PETSc/DM.pyx:1867

restoreGlobalVec(vg)#

Restore a global vector.

Not collective.

Parameters:

vg (Vec) – The global vector.

Return type:

None

Source code at petsc4py/PETSc/DM.pyx:809

restoreLocalVec(vl)#

Restore a local vector.

Not collective.

Parameters:

vl (Vec) – The local vector.

Return type:

None

Source code at petsc4py/PETSc/DM.pyx:842

setAppCtx(appctx)#

Set the application context.

Source code at petsc4py/PETSc/DM.pyx:339

Parameters:

appctx (Any)

Return type:

None

setAuxiliaryVec(aux, label, value=0, part=0)#

Set an auxiliary vector for a specific region.

Not collective.

Parameters:
  • aux (Vec) – The auxiliary vector.

  • label (DMLabel | None) – The name of the DMLabel.

  • value – Indicate the region.

  • part – The equation part, or 0 is unused.

Return type:

None

Source code at petsc4py/PETSc/DM.pyx:473

setBasicAdjacency(useCone, useClosure)#

Set the flags for determining variable influence.

Not collective.

Parameters:
  • useCone (bool) – Whether adjacency uses cone information.

  • useClosure (bool) – Whether adjacency is computed using full closure information.

Return type:

None

Source code at petsc4py/PETSc/DM.pyx:349

setCellCoordinateDM(dm)#

Set the cell coordinate DM.

Collective.

Parameters:

dm (DM) – The cell coordinate DM.

Return type:

None

Source code at petsc4py/PETSc/DM.pyx:1068

setCellCoordinateSection(dim, sec)#

Set the cell coordinate layout over the DM.

Collective.

Parameters:
Return type:

None

Source code at petsc4py/PETSc/DM.pyx:1100

setCellCoordinates(c)#

Set a global vector with the cellwise coordinates.

Collective.

Parameters:

c (Vec) – The global cell coordinate vector.

Return type:

None

Source code at petsc4py/PETSc/DM.pyx:1135

setCellCoordinatesLocal(c)#

Set a local vector with the cellwise coordinates.

Not collective.

Parameters:

c (Vec) – The local cell coordinate vector.

Return type:

None

Source code at petsc4py/PETSc/DM.pyx:1167

setCoarseDM(dm)#

Set the coarse DM.

Collective.

See also

DMSetCoarseDM

Source code at petsc4py/PETSc/DM.pyx:961

Parameters:

dm (DM)

Return type:

None

setCoordinateDim(dim)#

Set the dimension of embedding space for coordinates values.

Not collective.

Parameters:

dim (int) – The embedding dimension.

Return type:

None

Source code at petsc4py/PETSc/DM.pyx:252

setCoordinateDisc(disc, project)#

Project coordinates to a different space.

Collective.

Parameters:
  • disc (FE) – The new coordinates discretization.

  • project (bool)

Return type:

Self

Source code at petsc4py/PETSc/DM.pyx:1199

setCoordinates(c)#

Set a global vector that holds the coordinates.

Collective.

Parameters:

c (Vec) – Coordinate Vector.

Return type:

None

See also

DMSetCoordinates

Source code at petsc4py/PETSc/DM.pyx:1004

setCoordinatesLocal(c)#

Set a local vector with the ghost point holding the coordinates.

Not collective.

Parameters:

c (Vec) – Coordinate Vector.

Return type:

None

Source code at petsc4py/PETSc/DM.pyx:1036

setDimension(dim)#

Set the topological dimension of the DM.

Collective.

Parameters:

dim (int) – Topological dimension.

Return type:

None

See also

DMSetDimension

Source code at petsc4py/PETSc/DM.pyx:220

setField(index, field, label=None)#

Set the discretization object for a given DM field.

Logically collective.

Parameters:
  • index (int) – The field number.

  • field (Object) – The discretization object.

  • label (str | None) – The name of the label indicating the support of the field, or None for the entire mesh.

Return type:

None

See also

DMSetField

Source code at petsc4py/PETSc/DM.pyx:560

setFieldAdjacency(field, useCone, useClosure)#

Set the flags for determining variable influence.

Not collective.

Parameters:
  • field (int) – The field number.

  • useCone (bool) – Whether adjacency uses cone information.

  • useClosure (bool) – Whether adjacency is computed using full closure information.

Return type:

None

See also

DMSetAdjacency

Source code at petsc4py/PETSc/DM.pyx:392

setFromOptions()#

Configure the object from the options database.

Collective.

Source code at petsc4py/PETSc/DM.pyx:312

Return type:

None

setGlobalSection(sec)#

Set the Section encoding the global data layout for the DM.

Collective.

Source code at petsc4py/PETSc/DM.pyx:1676

Parameters:

sec (Section)

Return type:

None

setKSPComputeOperators(operators, args=None, kargs=None)#

Matrix associated with the linear system.

Collective.

Parameters:
  • operator – Callback function to compute the operators.

  • args (tuple[Any, ...] | None) – Positional arguments for the callback.

  • kargs (dict[str, Any] | None) – Keyword arguments for the callback.

Return type:

None

Source code at petsc4py/PETSc/DM.pyx:2120

setLabelOutput(name, output)#

Set if a given label should be saved to a view.

Not collective.

Parameters:
  • name (str) – The label name.

  • output (bool) – If True, the label is saved to the viewer.

Return type:

None

See also

DMSetLabelOutput

Source code at petsc4py/PETSc/DM.pyx:2072

setLabelValue(name, point, value)#

Set a point to a DMLabel with a given value.

Not collective.

Parameters:
  • name (str) – The label name.

  • point (int) – The mesh point.

  • value (int) – The label value for the point.

Return type:

None

See also

DMSetLabelValue

Source code at petsc4py/PETSc/DM.pyx:1912

setLocalSection(sec)#

Set the Section encoding the local data layout for the DM.

Collective.

Source code at petsc4py/PETSc/DM.pyx:1649

Parameters:

sec (Section)

Return type:

None

setMatType(mat_type)#

Set matrix type to be used by DM.createMat.

Logically collective.

Parameters:

mat_type (Type | str) – The matrix type.

Return type:

None

Notes

The option -dm_mat_type is used to set the matrix type.

Source code at petsc4py/PETSc/DM.pyx:1283

setNumFields(numFields)#

Set the number of fields in the DM.

Logically collective.

See also

DMSetNumFields

Source code at petsc4py/PETSc/DM.pyx:533

Parameters:

numFields (int)

Return type:

None

setOptionsPrefix(prefix)#

Set the prefix used for searching for options in the database.

Logically collective.

Source code at petsc4py/PETSc/DM.pyx:270

Parameters:

prefix (str | None)

Return type:

None

setPointSF(sf)#

Set the SF encoding the parallel DOF overlap for the DM.

Logically collective.

See also

DMSetPointSF

Source code at petsc4py/PETSc/DM.pyx:1781

Parameters:

sf (SF)

Return type:

None

setRefineLevel(level)#

Set the number of refinements.

Not collective.

Parameters:
  • nlevels – The number of refinement.

  • level (int)

Return type:

None

See also

DMSetRefineLevel

Source code at petsc4py/PETSc/DM.pyx:1537

setSNESFunction(function, args=None, kargs=None)#

Set SNES residual evaluation function.

Not collective.

Parameters:
  • function (SNESFunction) – The callback.

  • args (tuple[Any, ...] | None) – Positional arguments for the callback.

  • kargs (dict[str, Any] | None) – Keyword arguments for the callback.

Return type:

None

Source code at petsc4py/PETSc/DM.pyx:2196

setSNESJacobian(jacobian, args=None, kargs=None)#

Set the SNES Jacobian evaluation function.

Not collective.

Parameters:
  • jacobian (SNESJacobianFunction) – The Jacobian callback.

  • args (tuple[Any, ...] | None) – Positional arguments for the callback.

  • kargs (dict[str, Any] | None) – Keyword arguments for the callback.

Return type:

None

Source code at petsc4py/PETSc/DM.pyx:2228

setSectionSF(sf)#

Set the Section encoding the parallel DOF overlap for the DM.

Logically collective.

See also

DMSetSectionSF

Source code at petsc4py/PETSc/DM.pyx:1750

Parameters:

sf (SF)

Return type:

None

setType(dm_type)#

Build a DM.

Collective.

Parameters:

dm_type (Type | str) – The type of DM.

Return type:

None

Notes

DM types are available in DM.Type class.

See also

DM.Type, DMSetType

Source code at petsc4py/PETSc/DM.pyx:169

setUp()#

Return the data structure.

Collective.

See also

DMSetUp

Source code at petsc4py/PETSc/DM.pyx:324

Return type:

Self

setVecType(vec_type)#

Set the type of vector.

Logically collective.

Source code at petsc4py/PETSc/DM.pyx:752

Parameters:

vec_type (Type | str)

Return type:

None

view(viewer=None)#

View the DM.

Collective.

Parameters:

viewer (Viewer | None) – The DM viewer.

Return type:

None

See also

DMView

Source code at petsc4py/PETSc/DM.pyx:76

Attributes Documentation

appctx#

Application context.

Source code at petsc4py/PETSc/DM.pyx:2312

ds#

Discrete space.

Source code at petsc4py/PETSc/DM.pyx:2322