petsc4py.PETSc.DM#
- class petsc4py.PETSc.DM#
Bases:
Object
An object describing a computational grid or mesh.
Enumerations
DM
Boundary types.The
DM
cell types.The
DM
reordering default flags.DM
types.Methods Summary
adaptLabel
(label)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
.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)convert
(dm_type)copyDS
(dm[, minDegree, maxDegree])copyDisc
(dm)copyFields
(dm[, minDegree, maxDegree])create
([comm])Return an empty
DM
.createDS
()Create discrete systems.
Return a list of
IS
objects.Return a global vector.
createInjection
(dm)Return the injection matrix into a finer
DM
.Return the interpolation matrix to a finer
DM
.createLabel
(name)Create a label of the given name if it does not already exist.
Return a local vector.
createMassMatrix
(dmf)Return an empty matrix.
Return the restriction matrix between this
DM
and the givenDM
.createSectionSF
(localsec, globalsec)createSubDM
(fields)destroy
()Destroy the object.
Return the application context.
getAuxiliaryVec
([label, value, part])Return an auxiliary vector for region.
Return the flags for determining variable influence.
Return the inherent block size associated with a
DM
.Return the dimension of embedding space for coordinates values.
Return the cell coordinate
DM
.Return the cell coordinate layout over the
DM
.Return a global vector with the cellwise coordinates.
Return a local vector with the cellwise coordinates.
Return the coarse
DM
.Return the number of coarsenings.
Return the coordinate
DM
.Return the dimension of embedding space for coordinates values.
Return coordinate values layout over the mesh.
Return a global vector with the coordinates associated.
Return a local vector with the coordinates associated.
Check if the coordinates have been localized for cells.
getDS
()Return default
DS
.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.
Return the
Section
encoding the global data layout for theDM
.getGlobalVec
([name])Return a global vector.
getLGMap
()Return local mapping to global mapping.
getLabel
(name)Return the label of a given name.
getLabelIdIS
(name)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.Return the bounding box for the piece of the
DM
.Return the
Section
encoding the local data layout for theDM
.getLocalVec
([name])Return a local vector.
Return the number of fields in the
DM
.Return the number of labels defined by on the
DM
.Return the prefix used for searching for options in the database.
Return the refinement level.
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.
Create local coordinates for cells having periodic faces.
refine
([comm])Return a refined
DM
object.refineHierarchy
(nlevels)removeLabel
(name)Remove and destroy the label by name.
restoreGlobalVec
(vg[, name])Restore a global vector obtained with
getGlobalVec
.restoreLocalVec
(vl[, name])Restore a local vector obtained with
getLocalVec
.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.
Set the cell coordinate
DM
.setCellCoordinateSection
(dim, sec)Set the cell coordinate layout over the
DM
.Set a global vector with the cellwise coordinates.
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.
Set a global vector that holds the coordinates.
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.
Configure the object from the options database.
setGlobalSection
(sec)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)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)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 theDM
.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
Application context.
Discrete space.
Methods Documentation
- adaptMetric(metric, bdLabel=None, rgLabel=None)#
Return a mesh adapted to the specified metric field.
Collective.
- Parameters:
- Return type:
See also
- addCoarsenHook(coarsenhook, restricthook, args=None, kargs=None)#
Add a callback to be executed when restricting to a coarser grid.
Logically collective.
- Parameters:
coarsenhook (DMCoarsenHookFunction) – The coarsen hook function.
restricthook (DMRestrictHookFunction) – The restrict hook function.
args (tuple[Any, ...] | None) – Positional arguments for the hooks.
kargs (dict[str, Any] | None) – Keyword arguments for the hooks.
- Return type:
See also
- addField(field, label=None)#
Add a field to a
DM
object.Logically collective.
- Parameters:
- Return type:
See also
- appendOptionsPrefix(prefix)#
Append to the prefix used for searching for options in the database.
Logically collective.
- clearDS()#
Remove all discrete systems from the
DM
.Logically collective.
See also
Source code at petsc4py/PETSc/DM.pyx:691
- Return type:
- clearFields()#
Remove all fields from the
DM
.Logically collective.
See also
Source code at petsc4py/PETSc/DM.pyx:634
- Return type:
- clearLabelStratum(name, value)#
Remove all points from a stratum.
Not collective.
See also
- clearLabelValue(name, point, value)#
Remove a point from a
DMLabel
with given value.Not collective.
- Parameters:
- Return type:
See also
- clone()#
Return the cloned
DM
.Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:155
- Return type:
- coarsen(comm=None)#
Return a coarsened
DM
object.Collective.
- Parameters:
comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.- Return type:
See also
- coarsenHierarchy(nlevels)#
Coarsen this
DM
and return the coarsenedDM
hierarchy.Collective.
See also
- copyDS(dm, minDegree=None, maxDegree=None)#
Copy the discrete systems for this
DM
into anotherDM
.Collective.
- Parameters:
- Return type:
See also
- copyDisc(dm)#
Copy fields and discrete systems of a
DM
into anotherDM
.Collective.
See also
- copyFields(dm, minDegree=None, maxDegree=None)#
Copy the discretizations of this
DM
into anotherDM
.Collective.
- Parameters:
- Return type:
See also
- create(comm=None)#
Return an empty
DM
.Collective.
- Parameters:
comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.- Return type:
See also
- createDS()#
Create discrete systems.
Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:679
- Return type:
- createFieldDecomposition()#
Return a list of
IS
objects.Not collective.
Notes
The user is responsible for freeing all requested arrays.
See also
- createGlobalVec()#
Return a global vector.
Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:798
- Return type:
- createLabel(name)#
Create a label of the given name if it does not already exist.
Not collective.
See also
- createLocalVec()#
Return a local vector.
Not collective.
See also
Source code at petsc4py/PETSc/DM.pyx:812
- Return type:
- createMassMatrix(dmf)#
Return the mass matrix between this
DM
and the givenDM
.Collective.
See also
- createMat()#
Return an empty matrix.
Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1380
- Return type:
- createRestriction(dm)#
Return the restriction matrix between this
DM
and the givenDM
.Collective.
See also
- createSectionSF(localsec, globalsec)#
Create the
SF
encoding the parallel DOF overlap for theDM
.Collective.
- Parameters:
- Return type:
Notes
Encoding based on the
Section
describing the data layout.See also
- createSubDM(fields)#
Return
IS
andDM
encapsulating a subproblem.Not collective.
- Returns:
- Parameters:
- Return type:
See also
- destroy()#
Destroy the object.
Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:121
- Return type:
- getAppCtx()#
Return the application context.
Source code at petsc4py/PETSc/DM.pyx:343
- Return type:
- getAuxiliaryVec(label=None, value=0, part=0)#
Return an auxiliary vector for region.
Not collective.
- Parameters:
- Return type:
See also
- getBasicAdjacency()#
Return the flags for determining variable influence.
Not collective.
- Returns:
- Return type:
See also
- getBlockSize()#
Return the inherent block size associated with a
DM
.Not collective.
See also
Source code at petsc4py/PETSc/DM.pyx:770
- Return type:
- getBoundingBox()#
Return the dimension of embedding space for coordinates values.
Not collective.
See also
- getCellCoordinateDM()#
Return the cell coordinate
DM
.Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1159
- Return type:
- getCellCoordinateSection()#
Return the cell coordinate layout over the
DM
.Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1194
- Return type:
- getCellCoordinates()#
Return a global vector with the cellwise coordinates.
Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1226
- Return type:
- getCellCoordinatesLocal()#
Return a local vector with the cellwise coordinates.
Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1258
- Return type:
- getCoarseDM()#
Return the coarse
DM
.Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1020
- Return type:
- getCoarsenLevel()#
Return the number of coarsenings.
Not collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1629
- Return type:
- getCoordinateDM()#
Return the coordinate
DM
.Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1048
- Return type:
- getCoordinateDim()#
Return the dimension of embedding space for coordinates values.
Not collective.
See also
Source code at petsc4py/PETSc/DM.pyx:238
- Return type:
- getCoordinateSection()#
Return coordinate values layout over the mesh.
Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1063
- Return type:
- getCoordinates()#
Return a global vector with the coordinates associated.
Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1095
- Return type:
- getCoordinatesLocal()#
Return a local vector with the coordinates associated.
Collective the first time it is called.
See also
Source code at petsc4py/PETSc/DM.pyx:1127
- Return type:
- getCoordinatesLocalized()#
Check if the coordinates have been localized for cells.
Not collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1292
- Return type:
- getDS()#
Return default
DS
.Not collective.
See also
Source code at petsc4py/PETSc/DM.pyx:703
- Return type:
- getDimension()#
Return the topological dimension of the
DM
.Not collective.
See also
Source code at petsc4py/PETSc/DM.pyx:206
- Return type:
- getFieldAdjacency(field)#
Return the flags for determining variable influence.
Not collective.
- Parameters:
field (int) – The field number.
- Returns:
- Return type:
See also
- getGlobalSection()#
Return the
Section
encoding the global data layout for theDM
.Collective the first time it is called.
See also
Source code at petsc4py/PETSc/DM.pyx:1762
- Return type:
- getGlobalVec(name=None)#
Return a global vector.
Collective.
Notes
When done with the vector, it must be restored using
restoreGlobalVec
.
- getLGMap()#
Return local mapping to global mapping.
Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1003
- Return type:
- getLabel(name)#
Return the label of a given name.
Not collective.
See also
- getLabelName(index)#
Return the name of nth label.
Not collective.
See also
- getLabelOutput(name)#
Return the output flag for a given label.
Not collective.
See also
- getLabelValue(name, point)#
Return the value in
DMLabel
for the given point.Not collective.
See also
- getLocalSection()#
Return the
Section
encoding the local data layout for theDM
.Not collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1735
- Return type:
- getLocalVec(name=None)#
Return a local vector.
Not collective.
Notes
When done with the vector, it must be restored using
restoreLocalVec
.See also
- getNumFields()#
Return the number of fields in the
DM
.Not collective.
See also
Source code at petsc4py/PETSc/DM.pyx:546
- Return type:
- getNumLabels()#
Return the number of labels defined by on the
DM
.Not collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1867
- Return type:
- getOptionsPrefix()#
Return the prefix used for searching for options in the database.
Not collective.
Source code at petsc4py/PETSc/DM.pyx:284
- Return type:
- getPointSF()#
Return the
SF
encoding the parallel DOF overlap for theDM
.Not collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1840
- Return type:
- getRefineLevel()#
Return the refinement level.
Not collective.
See also
Source code at petsc4py/PETSc/DM.pyx:1597
- Return type:
- getSectionSF()#
Return the
Section
encoding the parallel DOF overlap.Collective the first time it is called.
See also
Source code at petsc4py/PETSc/DM.pyx:1809
- Return type:
- getStratumIS(name, value)#
Return the points in a label stratum.
Not collective.
See also
- getStratumSize(name, value)#
Return the number of points in a label stratum.
Not collective.
See also
- getType()#
Return the
DM
type name.Not collective.
See also
Source code at petsc4py/PETSc/DM.pyx:192
- Return type:
- globalToLocal(vg, vl, addv=None)#
Update local vectors from global vector.
Neighborwise collective.
- Parameters:
vg (Vec) – The global vector.
vl (Vec) – The local vector.
addv (InsertModeSpec | None) – Insertion mode.
- Return type:
See also
- load(viewer)#
Return a
DM
stored in binary.Collective.
- Parameters:
viewer (Viewer) – Viewer used to store the
DM
, likeViewer.Type.BINARY
orViewer.Type.HDF5
.- Return type:
Notes
When using
Viewer.Type.HDF5
format, one can save multipleDMPlex
meshes in a single HDF5 files. This in turn requires one to name theDMPlex
object withObject.setName
before saving it withDM.view
and before loading it withDM.load
for identification of the mesh object.See also
- localToGlobal(vl, vg, addv=None)#
Update global vectors from local vector.
Neighborwise collective.
- Parameters:
vl (Vec) – The local vector.
vg (Vec) – The global vector.
addv (InsertModeSpec | None) – Insertion mode.
- Return type:
See also
- localToLocal(vl, vlg, addv=None)#
Map the values from a local vector to another local vector.
Neighborwise collective.
- Parameters:
vl (Vec) – The local vector.
vlg (Vec) – The global vector.
addv (InsertModeSpec | None) – Insertion mode.
- Return type:
See also
- localizeCoordinates()#
Create local coordinates for cells having periodic faces.
Collective.
Notes
Used if the mesh is periodic.
See also
Source code at petsc4py/PETSc/DM.pyx:1340
- Return type:
- refine(comm=None)#
Return a refined
DM
object.Collective.
- Parameters:
comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.- Return type:
See also
- removeLabel(name)#
Remove and destroy the label by name.
Not collective.
See also
- restoreGlobalVec(vg, name=None)#
Restore a global vector obtained with
getGlobalVec
.Logically collective.
- Parameters:
- Return type:
- restoreLocalVec(vl, name=None)#
Restore a local vector obtained with
getLocalVec
.Not collective.
- Parameters:
- Return type:
- setAppCtx(appctx)#
Set the application context.
- setAuxiliaryVec(aux, label, value=0, part=0)#
Set an auxiliary vector for a specific region.
Not collective.
- Parameters:
- Return type:
See also
- setBasicAdjacency(useCone, useClosure)#
Set the flags for determining variable influence.
Not collective.
- Parameters:
- Return type:
See also
- setCellCoordinateSection(dim, sec)#
Set the cell coordinate layout over the
DM
.Collective.
- Parameters:
- Return type:
See also
- setCellCoordinates(c)#
Set a global vector with the cellwise coordinates.
Collective.
See also
- setCellCoordinatesLocal(c)#
Set a local vector with the cellwise coordinates.
Not collective.
See also
- setCoordinateDim(dim)#
Set the dimension of embedding space for coordinates values.
Not collective.
See also
- setCoordinateDisc(disc, project)#
Project coordinates to a different space.
Collective.
See also
- setCoordinates(c)#
Set a global vector that holds the coordinates.
Collective.
See also
- setCoordinatesLocal(c)#
Set a local vector with the ghost point holding the coordinates.
Not collective.
See also
- setField(index, field, label=None)#
Set the discretization object for a given
DM
field.Logically collective.
- Parameters:
- Return type:
See also
- setFieldAdjacency(field, useCone, useClosure)#
Set the flags for determining variable influence.
Not collective.
- Parameters:
- Return type:
See also
- setFromOptions()#
Configure the object from the options database.
Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:312
- Return type:
- setGlobalSection(sec)#
Set the
Section
encoding the global data layout for theDM
.Collective.
See also
- setKSPComputeOperators(operators, args=None, kargs=None)#
Matrix associated with the linear system.
Collective.
- Parameters:
- Return type:
See also
- setLabelOutput(name, output)#
Set if a given label should be saved to a view.
Not collective.
- Parameters:
- Return type:
See also
- setLabelValue(name, point, value)#
Set a point to a
DMLabel
with a given value.Not collective.
- Parameters:
- Return type:
See also
- setLocalSection(sec)#
Set the
Section
encoding the local data layout for theDM
.Collective.
See also
- setMatType(mat_type)#
Set matrix type to be used by
DM.createMat
.Logically collective.
Notes
The option
-dm_mat_type
is used to set the matrix type.See also
- setOptionsPrefix(prefix)#
Set the prefix used for searching for options in the database.
Logically collective.
- setPointSF(sf)#
Set the
SF
encoding the parallel DOF overlap for theDM
.Logically collective.
See also
- setRefineLevel(level)#
Set the number of refinements.
Not collective.
See also
- setSNESFunction(function, args=None, kargs=None)#
Set
SNES
residual evaluation function.Not collective.
- Parameters:
- Return type:
See also
- setSNESJacobian(jacobian, args=None, kargs=None)#
Set the
SNES
Jacobian evaluation function.Not collective.
- Parameters:
- Return type:
See also
- setSectionSF(sf)#
Set the
Section
encoding the parallel DOF overlap for theDM
.Logically collective.
See also
- setUp()#
Return the data structure.
Collective.
See also
Source code at petsc4py/PETSc/DM.pyx:324
- Return type:
- setVecType(vec_type)#
Set the type of vector.
Logically collective.
See also
Attributes Documentation
- appctx#
Application context.
- ds#
Discrete space.