petsc4py.PETSc.DMPlex#
- class petsc4py.PETSc.DMPlex#
Bases:
DM
Encapsulate an unstructured mesh.
DMPlex encapsulates both topology and geometry. It is capable of parallel refinement and coarsening (using Pragmatic or ParMmg) and parallel redistribution for load balancing.
Methods Summary
computeCellGeometryFVM
(cell)Compute the volume for a given cell.
computeGradientClementInterpolant
(locX, locC)Return the L2 projection of the cellwise gradient of a function onto P1.
constructGhostCells
([labelName])Construct ghost cells which connect to every boundary face.
coordinatesLoad
(viewer, sfxc)Load coordinates into this
DMPlex
object.coordinatesView
(viewer)Save
DMPlex
coordinates into a file.create
([comm])Create a
DMPlex
object, which encapsulates an unstructured mesh.createBoxMesh
(faces[, lower, upper, ...])Create a mesh on the tensor product of intervals.
createBoxSurfaceMesh
(faces[, lower, upper, ...])Create a mesh on the surface of a box mesh using tensor cells.
createCGNS
(cgid[, interpolate, comm])Create a
DMPlex
mesh from a CGNS file.createCGNSFromFile
(filename[, interpolate, comm])"Create a
DMPlex
mesh from a CGNS file.createClosureIndex
(sec)Calculate an index for
sec
for the closure operation.Create an
IS
covering the coarseDMPlex
chart with the fine points as data.createCohesiveSubmesh
(hasLagrange, value)Extract the hypersurface defined by one face of the cohesive cells.
createExodus
(exoid[, interpolate, comm])Create a
DMPlex
mesh from an ExodusII file ID.createExodusFromFile
(filename[, ...])Create a
DMPlex
mesh from an ExodusII file.createFromCellList
(dim, cells, coords[, ...])Create a
DMPlex
from a list of vertices for each cell on process 0.createFromFile
(filename[, plexname, ...])Create
DMPlex
from a file.createGmsh
(viewer[, interpolate, comm])Create a
DMPlex
mesh from a Gmsh file viewer.Create a global numbering for all points.
createSection
(numComp, numDof[, bcField, ...])Create a
Section
based upon the DOF layout specification provided.distribute
([overlap])Distribute the mesh and any associated sections.
distributeField
(sf, sec, vec[, newsec, newvec])Distribute field data with a with a given
SF
.Return a flag indicating whether the
DM
should be distributed by default.distributeOverlap
([overlap])Add partition overlap to a distributed non-overlapping
DMPlex
.distributeSetDefault
(flag)Set flag indicating whether the
DMPlex
should be distributed by default.Retrieve the name of the specific parallel distribution.
distributionSetName
(name)Set the name of the specific parallel distribution.
generate
(boundary[, name, interpolate])Generate a mesh.
getAdjacency
(p)Return all points adjacent to the given point.
Query whether adjacency in the mesh uses the point-to-point constraints.
Return a global cell numbering for all cells on this process.
getCellType
(p)Return the polytope type of a given cell.
Return the
DMLabel
recording the polytope type of each cell.getChart
()Return the interval for all mesh points [
pStart
,pEnd
).getCone
(p)Return the points on the in-edges for this point in the DAG.
Return the orientations on the in-edges for this point in the DAG.
getConeSize
(p)Return the number of in-edges for this point in the DAG.
getDepth
()Return the depth of the DAG representing this mesh.
getDepthStratum
(svalue)Return the bounds [
start
,end
) for all points at a certain depth.getFullJoin
(points)Return an array for the join of the set of points.
getHeightStratum
(svalue)Return the bounds [
start
,end
) for all points at a certain height.getJoin
(points)Return an array for the join of the set of points.
Return the maximum number of in-edges and out-edges of the DAG.
getMeet
(points)Return an array for the meet of the set of points.
Return the minimum distance from any cell centroid to a face.
getOrdering
(otype)Calculate a reordering of the mesh.
Return the mesh partitioner.
getPointDepth
(point)Return the depth of a given point.
getPointGlobal
(point)Return location of point data in global
Vec
.getPointGlobalField
(point, field)Return location of point field data in global
Vec
.getPointHeight
(point)Return the height of a given point.
getPointLocal
(point)Return location of point data in local
Vec
.getPointLocalField
(point, field)Return location of point field data in local
Vec
.Retrieve the maximum cell volume for refinement.
Retrieve the flag for uniform refinement.
Return an
IS
covering the entire subdm chart.Return a
DMLabel
with point dimension as values.getSupport
(p)Return the points on the out-edges for this point in the DAG.
Return the number of out-edges for this point in the DAG.
getTransitiveClosure
(p[, useCone])Return the points and orientations on the transitive closure of this point.
getVecClosure
(sec, vec, point)Return an array of the values on the closure of a point.
Return a global vertex numbering for all vertices on this process.
globalVectorLoad
(viewer, sectiondm, sf, vec)Load on-disk vector data into a global vector.
globalVectorView
(viewer, sectiondm, vec)Save a global vector.
insertCone
(p, conePos, conePoint)DMPlexInsertCone - Insert a point into the in-edges for the point p in the DAG.
insertConeOrientation
(p, conePos, ...)Insert a point orientation for the in-edge for the point p in the DAG.
Convert to a mesh with all intermediate faces, edges, etc.
Return the flag indicating if the mesh is distributed.
Return the flag indicating if the first cell is a simplex.
labelCohesiveComplete
(label, bdlabel, ...)Add all other mesh pieces to complete the surface.
labelComplete
(label)Add the transitive closure to the surface.
labelsLoad
(viewer, sfxc)Load labels into this
DMPlex
object.labelsView
(viewer)Save
DMPlex
labels into a file.localVectorLoad
(viewer, sectiondm, sf, vec)Load on-disk vector data into a local vector.
localVectorView
(viewer, sectiondm, vec)Save a local vector.
markBoundaryFaces
(label[, value])Mark all faces on the boundary.
metricAverage2
(metric1, metric2, metricAvg)Compute and return the unweighted average of two metrics.
metricAverage3
(metric1, metric2, metric3, ...)Compute and return the unweighted average of three metrics.
metricCreate
([field])Create a Riemannian metric field.
metricCreateIsotropic
(indicator[, field])Construct an isotropic metric from an error indicator.
metricCreateUniform
(alpha[, field])Construct a uniform isotropic metric.
metricDeterminantCreate
([field])Create the determinant field for a Riemannian metric.
metricEnforceSPD
(metric, ometric, determinant)Enforce symmetric positive-definiteness of a metric.
Return the metric gradation factor.
Return the metric Hausdorff number.
Return the maximum tolerated metric anisotropy.
Return the maximum tolerated metric magnitude.
Return the minimum tolerated metric magnitude.
Return the order p for L-p normalization.
Return the number of parallel adaptation iterations.
Return the target metric complexity.
Return the verbosity of the mesh adaptation package.
metricIntersection2
(metric1, metric2, metricInt)Compute and return the intersection of two metrics.
metricIntersection3
(metric1, metric2, ...)Compute the intersection of three metrics.
Return the flag indicating whether the metric is isotropic or not.
Return the flag indicating whether the metric is uniform or not.
Return the flag indicating whether node insertion and deletion are turned off.
Return the flag indicating whether node movement is turned off.
Return the flag indicating whether surface modification is turned off.
Return the flag indicating whether facet swapping is turned off.
metricNormalize
(metric, ometric, determinant)Apply L-p normalization to a metric.
Return
true
if anisotropy is restricted before normalization.Configure the object from the options database.
metricSetGradationFactor
(beta)Set the metric gradation factor.
metricSetHausdorffNumber
(hausd)Set the metric Hausdorff number.
metricSetIsotropic
(isotropic)Record whether the metric is isotropic or not.
metricSetMaximumAnisotropy
(a_max)Set the maximum tolerated metric anisotropy.
metricSetMaximumMagnitude
(h_max)Set the maximum tolerated metric magnitude.
metricSetMinimumMagnitude
(h_min)Set the minimum tolerated metric magnitude.
metricSetNoInsertion
(noInsert)Set the flag indicating whether node insertion should be turned off.
metricSetNoMovement
(noMove)Set the flag indicating whether node movement should be turned off.
metricSetNoSurf
(noSurf)Set the flag indicating whether surface modification should be turned off.
metricSetNoSwapping
(noSwap)Set the flag indicating whether facet swapping should be turned off.
Set the order p for L-p normalization.
metricSetNumIterations
(numIter)Set the number of parallel adaptation iterations.
Record whether anisotropy is be restricted before normalization or after.
metricSetTargetComplexity
(targetComplexity)Set the target metric complexity.
metricSetUniform
(uniform)Record whether the metric is uniform or not.
metricSetVerbosity
(verbosity)Set the verbosity of the mesh adaptation package.
orient
()Give a consistent orientation to the input mesh.
permute
(perm)Reorder the mesh according to the input permutation.
rebalanceSharedPoints
([entityDepth, ...])Redistribute shared points in order to achieve better balancing.
Return flag indicating whether the
DMPlex
should be reordered by default.reorderSetDefault
(flag)Set flag indicating whether the DM should be reordered by default.
sectionLoad
(viewer, sectiondm, sfxc)Load section into a
DM
.sectionView
(viewer, sectiondm)Save a section associated with a
DMPlex
.setAdjacencyUseAnchors
([useAnchors])Define adjacency in the mesh using the point-to-point constraints.
setCellType
(p, ctype)Set the polytope type of a given cell.
setChart
(pStart, pEnd)Set the interval for all mesh points [
pStart
,pEnd
).setCone
(p, cone[, orientation])Set the points on the in-edges for this point in the DAG.
setConeOrientation
(p, orientation)Set the orientations on the in-edges for this point in the DAG.
setConeSize
(p, size)Set the number of in-edges for this point in the DAG.
setMatClosure
(sec, gsec, mat, point, values)Set an array of the values on the closure of
point
.setPartitioner
(part)Set the mesh partitioner.
setRefinementLimit
(refinementLimit)Set the maximum cell volume for refinement.
setRefinementUniform
([refinementUniform])Set the flag for uniform refinement.
setSupport
(p, supp)Set the points on the out-edges for this point in the DAG.
setSupportSize
(p, size)Set the number of out-edges for this point in the DAG.
setTetGenOptions
(opts)Set the options used for the Tetgen mesh generator.
setTriangleOptions
(opts)Set the options used for the Triangle mesh generator.
setVecClosure
(sec, vec, point, values[, addv])Set an array of the values on the closure of
point
.stratify
()Calculate the strata of DAG.
Create support (out-edge) information from cone (in-edge) information.
topologyLoad
(viewer)Load a topology into this
DMPlex
object.topologyView
(viewer)Save a
DMPlex
topology into a file.Convert to a mesh with only cells and vertices.
vecGetClosure
(sec, vec, p)Return an array of values on the closure of
p
.Methods Documentation
- computeCellGeometryFVM(cell)#
Compute the volume for a given cell.
Not collective.
- Parameters:
cell (int) – The cell.
- Returns:
- Return type:
- computeGradientClementInterpolant(locX, locC)#
Return the L2 projection of the cellwise gradient of a function onto P1.
Collective.
- Parameters:
- Return type:
See also
- constructGhostCells(labelName=None)#
Construct ghost cells which connect to every boundary face.
Collective.
- Parameters:
labelName (str | None) – The name of the label specifying the boundary faces. Defaults to
"Face Sets"
.- Returns:
numGhostCells – The number of ghost cells added to the
DMPlex
.- Return type:
See also
- coordinatesLoad(viewer, sfxc)#
Load coordinates into this
DMPlex
object.Collective.
- Parameters:
sfxc (SF) – The
SF
returned bytopologyLoad
.
- Return type:
See also
DM
,DMPlex
,DM.load
,DMPlex.topologyLoad
,DMPlex.labelsLoad
,DM.view
,SF
,Viewer
,DMPlexCoordinatesLoad
- create(comm=None)#
Create a
DMPlex
object, which encapsulates an unstructured mesh.Collective.
- Parameters:
comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.- Return type:
See also
- createBoxMesh(faces, lower=(0, 0, 0), upper=(1, 1, 1), simplex=True, periodic=False, interpolate=True, comm=None)#
Create a mesh on the tensor product of intervals.
Collective.
- Parameters:
faces (Sequence[int]) – Number of faces per dimension, or
None
for the default.simplex (bool | None) –
True
for simplices,False
for tensor cells.periodic (Sequence | str | int | bool | None) – The boundary type for the X, Y, Z direction, or
None
forDM.BoundaryType.NONE
.interpolate (bool | None) – Flag to create intermediate mesh entities (edges, faces).
comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.
- Return type:
See also
DM
,DMPlex
,DM.setFromOptions
,DMPlex.createFromFile
,DM.setType
,DM.create
,DMPlexCreateBoxMesh
- createBoxSurfaceMesh(faces, lower=(0, 0, 0), upper=(1, 1, 1), interpolate=True, comm=None)#
Create a mesh on the surface of a box mesh using tensor cells.
Collective.
- createCGNS(cgid, interpolate=True, comm=None)#
Create a
DMPlex
mesh from a CGNS file.Collective.
- Parameters:
- Return type:
- createCGNSFromFile(filename, interpolate=True, comm=None)#
“Create a
DMPlex
mesh from a CGNS file.Collective.
- Parameters:
- Return type:
See also
DM
,DMPlex
,DMPlex.create
,DMPlex.createCGNS
,DMPlex.createExodus
,DMPlexCreateCGNS
- createClosureIndex(sec)#
Calculate an index for
sec
for the closure operation.Not collective.
- Parameters:
sec (Section) – The
Section
describing the layout in the local vector, orNone
to use the default section.- Return type:
See also
DM
,DMPlex
,Section
,DMPlex.vecGetClosure
,DMPlexCreateClosureIndex
- createCoarsePointIS()#
Create an
IS
covering the coarseDMPlex
chart with the fine points as data.Collective.
- Returns:
fpointIS – The
IS
of all the fine points which exist in the original coarse mesh.- Return type:
See also
DM
,DMPlex
,IS
,DM.refine
,DMPlex.setRefinementUniform
,DMPlexCreateCoarsePointIS
- createCohesiveSubmesh(hasLagrange, value)#
Extract the hypersurface defined by one face of the cohesive cells.
Collective.
- Parameters:
- Return type:
See also
- createExodus(exoid, interpolate=True, comm=None)#
Create a
DMPlex
mesh from an ExodusII file ID.Collective.
- Parameters:
- Return type:
See also
- createExodusFromFile(filename, interpolate=True, comm=None)#
Create a
DMPlex
mesh from an ExodusII file.Collective.
- Parameters:
- Return type:
See also
DM
,DMPlex
,DM.create
,DMPlex.createExodus
,DMPlexCreateExodusFromFile
- createFromCellList(dim, cells, coords, interpolate=True, comm=None)#
Create a
DMPlex
from a list of vertices for each cell on process 0.Collective.
- Parameters:
dim (int) – The topological dimension of the mesh.
cells (Sequence[int]) – An array of number of cells times number of vertices on each cell.
coords (Sequence[float]) – An array of number of vertices times spatial dimension for coordinates.
comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.
- Return type:
See also
DM
,DMPlex
,DMPlex.create
,DMPlex.interpolate
,DMPlexCreateFromCellListPetsc
- createFromFile(filename, plexname='unnamed', interpolate=True, comm=None)#
Create
DMPlex
from a file.Collective.
- Parameters:
- Return type:
- createGmsh(viewer, interpolate=True, comm=None)#
Create a
DMPlex
mesh from a Gmsh file viewer.Collective.
- Parameters:
- Return type:
Notes
-dm_plex_gmsh_hybrid
forces triangular prisms to use tensor order.-dm_plex_gmsh_periodic
allows for reading Gmsh periodic section.-dm_plex_gmsh_highorder
allows for generating high-order coordinates.-dm_plex_gmsh_project
projects high-order coordinates to a different space, use the prefix-dm_plex_gmsh_project_
to define the space.-dm_plex_gmsh_use_regions
generates labels with region names.-dm_plex_gmsh_mark_vertices
adds vertices to generated labels.-dm_plex_gmsh_multiple_tags
allows multiple tags for default labels.-dm_plex_gmsh_spacedim <d>
embedding space dimension.See also
DM
,DMPlex
,DM.create
, Working with PETSc options,DMPlexCreateGmsh
- createPointNumbering()#
Create a global numbering for all points.
Collective.
See also
DM
,DMPlex
,DMPlex.getCellNumbering
,DMPlexCreatePointNumbering
Source code at petsc4py/PETSc/DMPlex.pyx:907
- Return type:
- createSection(numComp, numDof, bcField=None, bcComps=None, bcPoints=None, perm=None)#
Create a
Section
based upon the DOF layout specification provided.Not collective.
- Parameters:
numComp (Sequence[int]) – An array of size
numFields
holding the number of components per field.numDof (Sequence[int]) – An array of size
numFields*(dim+1)
holding the number of DOFs per field on a mesh piece of dimensiondim
.bcField (Sequence[int] | None) – An array of size
numBC
giving the field number for each boundary condition, wherenumBC
is the number of boundary conditions.bcComps (Sequence[IS] | None) – An array of size
numBC
giving anIS
holding the field components to which each boundary condition applies.bcPoints (Sequence[IS] | None) – An array of size
numBC
giving anIS
holding theDMPlex
points to which each boundary condition applies.
- Return type:
See also
DM
,DMPlex
,DMPlex.create
,Section.create
,Section.setPermutation
,DMPlexCreateSection
- distribute(overlap=0)#
Distribute the mesh and any associated sections.
Collective.
- Parameters:
- Returns:
sf – The
SF
used for point distribution, orNone
if not distributed.- Return type:
See also
- distributeField(sf, sec, vec, newsec=None, newvec=None)#
Distribute field data with a with a given
SF
.Collective.
- Parameters:
- Returns:
- Return type:
See also
- distributeGetDefault()#
Return a flag indicating whether the
DM
should be distributed by default.Not collective.
- Returns:
dist – Flag indicating whether the
DMPlex
should be distributed by default.- Return type:
- distributeOverlap(overlap=0)#
Add partition overlap to a distributed non-overlapping
DMPlex
.Collective.
- Parameters:
overlap (int | None) – The overlap of partitions (the same on all ranks).
- Returns:
sf – The
SF
used for point distribution.- Return type:
See also
DM
,DMPlex
,SF
,DMPlex.create
,DMPlex.distribute
,DMPlexDistributeOverlap
- distributeSetDefault(flag)#
Set flag indicating whether the
DMPlex
should be distributed by default.Logically collective.
- distributionGetName()#
Retrieve the name of the specific parallel distribution.
Not collective.
- Returns:
name – The name of the specific parallel distribution.
- Return type:
- distributionSetName(name)#
Set the name of the specific parallel distribution.
Logically collective.
- generate(boundary, name=None, interpolate=True)#
Generate a mesh.
Not collective.
- Parameters:
- Return type:
See also
DM
,DMPlex
,DMPlex.create
,DM.refine
, Working with PETSc options,DMPlexGenerate
- getAdjacency(p)#
Return all points adjacent to the given point.
Not collective.
See also
- getAdjacencyUseAnchors()#
Query whether adjacency in the mesh uses the point-to-point constraints.
Not collective.
Source code at petsc4py/PETSc/DMPlex.pyx:1439
- Return type:
- getCellNumbering()#
Return a global cell numbering for all cells on this process.
Collective the first time it is called.
See also
DM
,DMPlex
,DMPlex.getVertexNumbering
,DMPlexGetCellNumbering
Source code at petsc4py/PETSc/DMPlex.pyx:877
- Return type:
- getCellType(p)#
Return the polytope type of a given cell.
Not collective.
- Parameters:
p (int) – The cell.
- Return type:
- getCellTypeLabel()#
Return the
DMLabel
recording the polytope type of each cell.Not collective.
See also
DM
,DMPlex
,DMPlex.getCellType
,DM.createLabel
,DMPlexGetCellTypeLabel
Source code at petsc4py/PETSc/DMPlex.pyx:698
- Return type:
- getChart()#
Return the interval for all mesh points [
pStart
,pEnd
).Not collective.
- Returns:
- Return type:
See also
- getCone(p)#
Return the points on the in-edges for this point in the DAG.
Not collective.
- Parameters:
p (int) – The point, which must lie in the chart set with
DMPlex.setChart
.- Return type:
See also
DM
,DMPlex
,DMPlex.getConeSize
,DMPlex.setCone
,DMPlex.setChart
,DMPlexGetCone
- getConeOrientation(p)#
Return the orientations on the in-edges for this point in the DAG.
Not collective.
- Parameters:
p (int) – The point, which must lie in the chart set with
DMPlex.setChart
.- Return type:
- getConeSize(p)#
Return the number of in-edges for this point in the DAG.
Not collective.
- Parameters:
p (int) – The point, which must lie in the chart set with
DMPlex.setChart
.- Return type:
See also
DM
,DMPlex
,DMPlex.create
,DMPlex.setConeSize
,DMPlex.setChart
,DMPlexGetConeSize
- getDepth()#
Return the depth of the DAG representing this mesh.
Not collective.
See also
DM
,DMPlex
,DMPlex.getDepthStratum
,DMPlex.symmetrize
,DMPlexGetDepth
Source code at petsc4py/PETSc/DMPlex.pyx:921
- Return type:
- getDepthStratum(svalue)#
Return the bounds [
start
,end
) for all points at a certain depth.Not collective.
- Parameters:
svalue (int) – The requested depth.
- Returns:
- Return type:
- getFullJoin(points)#
Return an array for the join of the set of points.
Not collective.
See also
DM
,DMPlex
,DMPlex.getJoin
,DMPlex.getMeet
,DMPlexGetFullJoin
- getHeightStratum(svalue)#
Return the bounds [
start
,end
) for all points at a certain height.Not collective.
- Parameters:
svalue (int) – The requested height.
- Returns:
- Return type:
See also
DM
,DMPlex
,DMPlex.getDepthStratum
,DMPlex.getDepth
,DMPlexGetHeightStratum
- getJoin(points)#
Return an array for the join of the set of points.
Not collective.
See also
- getMaxSizes()#
Return the maximum number of in-edges and out-edges of the DAG.
Not collective.
- Returns:
- Return type:
See also
DM
,DMPlex
,DMPlex.create
,DMPlex.setConeSize
,DMPlex.setChart
,DMPlexGetMaxSizes
- getMeet(points)#
Return an array for the meet of the set of points.
Not collective.
See also
- getMinRadius()#
Return the minimum distance from any cell centroid to a face.
Not collective.
See also
Source code at petsc4py/PETSc/DMPlex.pyx:1790
- Return type:
- getOrdering(otype)#
Calculate a reordering of the mesh.
Collective.
- Parameters:
otype (OrderingType) – Type of reordering, see
Mat.OrderingType
.- Returns:
perm – The point permutation.
- Return type:
See also
DMPlex
,DMPlex.permute
,Mat.OrderingType
,Mat.getOrdering
,DMPlexGetOrdering
- getPartitioner()#
Return the mesh partitioner.
Not collective.
See also
DM
,DMPlex
,Partitioner
,Section
,DMPlex.distribute
,DMPlex.setPartitioner
,Partitioner.create
,PetscPartitionerDMPlexPartition
,DMPlexGetPartitioner
Source code at petsc4py/PETSc/DMPlex.pyx:1497
- Return type:
- getPointDepth(point)#
Return the depth of a given point.
Not collective.
See also
DM
,DMPlex
,DMPlex.getDepthStratum
,DMPlex.getDepth
,DMPlexGetPointDepth
- getPointGlobal(point)#
Return location of point data in global
Vec
.Not collective.
- Parameters:
point (int) – The topological point.
- Returns:
- Return type:
- getPointGlobalField(point, field)#
Return location of point field data in global
Vec
.Not collective.
- Parameters:
- Returns:
- Return type:
- getPointHeight(point)#
Return the height of a given point.
Not collective.
See also
- getPointLocal(point)#
Return location of point data in local
Vec
.Not collective.
- Parameters:
point (int) – The topological point.
- Returns:
- Return type:
- getPointLocalField(point, field)#
Return location of point field data in local
Vec
.Not collective.
- Parameters:
- Returns:
- Return type:
- getRefinementLimit()#
Retrieve the maximum cell volume for refinement.
Not collective.
See also
DM
,DMPlex
,DM.refine
,DMPlex.setRefinementLimit
,DMPlex.getRefinementUniform
,DMPlex.setRefinementUniform
,DMPlexGetRefinementLimit
Source code at petsc4py/PETSc/DMPlex.pyx:2101
- Return type:
- getRefinementUniform()#
Retrieve the flag for uniform refinement.
Not collective.
- Returns:
refinementUniform – The flag for uniform refinement.
- Return type:
- getSupport(p)#
Return the points on the out-edges for this point in the DAG.
Not collective.
- Parameters:
p (int) – The point, which must lie in the chart set with
DMPlex.setChart
.- Return type:
- getSupportSize(p)#
Return the number of out-edges for this point in the DAG.
Not collective.
- Parameters:
p (int) – The point, which must lie in the chart set with
DMPlex.setChart
.- Return type:
- getTransitiveClosure(p, useCone=True)#
Return the points and orientations on the transitive closure of this point.
Not collective.
- Parameters:
- Returns:
- Return type:
- getVecClosure(sec, vec, point)#
Return an array of the values on the closure of a point.
Not collective.
- Parameters:
- Return type:
See also
- getVertexNumbering()#
Return a global vertex numbering for all vertices on this process.
Collective the first time it is called.
See also
DM
,DMPlex
,DMPlex.getCellNumbering
,DMPlexGetVertexNumbering
Source code at petsc4py/PETSc/DMPlex.pyx:892
- Return type:
- globalVectorLoad(viewer, sectiondm, sf, vec)#
Load on-disk vector data into a global vector.
Collective.
- Parameters:
- Return type:
- globalVectorView(viewer, sectiondm, vec)#
Save a global vector.
Collective.
- Parameters:
- Return type:
- insertCone(p, conePos, conePoint)#
DMPlexInsertCone - Insert a point into the in-edges for the point p in the DAG.
Not collective.
- Parameters:
p (int) – The point, which must lie in the chart set with
DMPlex.setChart
.conePos (int) – The local index in the cone where the point should be put.
conePoint (int) – The mesh point to insert.
- Return type:
See also
DM
,DMPlex
,DMPlex.create
,DMPlex.getCone
,DMPlex.setChart
,DMPlex.setConeSize
,DM.setUp
,DMPlexInsertCone
- insertConeOrientation(p, conePos, coneOrientation)#
Insert a point orientation for the in-edge for the point p in the DAG.
Not collective.
- Parameters:
p (int) – The point, which must lie in the chart set with
DMPlex.setChart
conePos (int) – The local index in the cone where the point should be put.
coneOrientation (int) – The point orientation to insert.
- Return type:
- interpolate()#
Convert to a mesh with all intermediate faces, edges, etc.
Collective.
Source code at petsc4py/PETSc/DMPlex.pyx:1715
- Return type:
- isDistributed()#
Return the flag indicating if the mesh is distributed.
Collective.
See also
Source code at petsc4py/PETSc/DMPlex.pyx:1605
- Return type:
- isSimplex()#
Return the flag indicating if the first cell is a simplex.
Not collective.
See also
DM
,DMPlex
,DMPlex.getCellType
,DMPlex.getHeightStratum
,DMPlexIsSimplex
Source code at petsc4py/PETSc/DMPlex.pyx:1619
- Return type:
- labelCohesiveComplete(label, bdlabel, bdvalue, flip, subdm)#
Add all other mesh pieces to complete the surface.
Not collective.
- Parameters:
- Return type:
See also
DM
,DMPlex
,DMPlex.labelComplete
,DMPlexLabelCohesiveComplete
- labelComplete(label)#
Add the transitive closure to the surface.
Not collective.
See also
DM
,DMPlex
,DMPlex.labelCohesiveComplete
,DMPlexLabelComplete
- labelsLoad(viewer, sfxc)#
Load labels into this
DMPlex
object.Collective.
- Parameters:
sfxc (SF) – The
SF
returned bytopologyLoad
.
- Return type:
See also
DM
,DMPlex
,DM.load
,DMPlex.topologyLoad
,DMPlex.coordinatesLoad
,DM.view
,SF
,Viewer
,DMPlexLabelsLoad
- localVectorLoad(viewer, sectiondm, sf, vec)#
Load on-disk vector data into a local vector.
Collective.
- Parameters:
- Return type:
- localVectorView(viewer, sectiondm, vec)#
Save a local vector.
Collective.
- Parameters:
- Return type:
- markBoundaryFaces(label, value=None)#
Mark all faces on the boundary.
Not collective.
- Parameters:
- Return type:
See also
DM
,DMPlex
,DMLabel.create
,DM.createLabel
,DMPlexMarkBoundaryFaces
- metricAverage2(metric1, metric2, metricAvg)#
Compute and return the unweighted average of two metrics.
Collective.
- Parameters:
- Return type:
See also
- metricAverage3(metric1, metric2, metric3, metricAvg)#
Compute and return the unweighted average of three metrics.
Collective.
- Parameters:
- Return type:
See also
- metricCreate(field=0)#
Create a Riemannian metric field.
Collective.
- metricCreateIsotropic(indicator, field=0)#
Construct an isotropic metric from an error indicator.
Collective.
- Parameters:
- Return type:
- metricCreateUniform(alpha, field=0)#
Construct a uniform isotropic metric.
Collective.
- Parameters:
- Return type:
- metricDeterminantCreate(field=0)#
Create the determinant field for a Riemannian metric.
Collective.
- Parameters:
- Returns:
- Return type:
- metricEnforceSPD(metric, ometric, determinant, restrictSizes=False, restrictAnisotropy=False)#
Enforce symmetric positive-definiteness of a metric.
Collective.
- Parameters:
metric (Vec) – The metric.
ometric (Vec) – The output metric.
determinant (Vec) – The output determinant.
restrictSizes (bool | None) – Flag indicating whether maximum/minimum magnitudes should be enforced.
restrictAnisotropy (bool | None) – Flag indicating whether maximum anisotropy should be enforced.
- Returns:
- Return type:
- metricGetGradationFactor()#
Return the metric gradation factor.
Not collective.
See also
DMPlex.metricSetGradationFactor
,DMPlex.metricGetHausdorffNumber
,DMPlexMetricGetGradationFactor
Source code at petsc4py/PETSc/DMPlex.pyx:2828
- Return type:
- metricGetHausdorffNumber()#
Return the metric Hausdorff number.
Not collective.
See also
DMPlex.metricGetGradationFactor
,DMPlex.metricSetHausdorffNumber
,DMPlexMetricGetHausdorffNumber
Source code at petsc4py/PETSc/DMPlex.pyx:2862
- Return type:
- metricGetMaximumAnisotropy()#
Return the maximum tolerated metric anisotropy.
Not collective.
See also
DMPlex.metricSetMaximumAnisotropy
,DMPlex.metricGetMaximumMagnitude
,DMPlexMetricGetMaximumAnisotropy
Source code at petsc4py/PETSc/DMPlex.pyx:2726
- Return type:
- metricGetMaximumMagnitude()#
Return the maximum tolerated metric magnitude.
Not collective.
See also
DMPlex.metricSetMaximumMagnitude
,DMPlex.metricGetMinimumMagnitude
,DMPlexMetricGetMaximumMagnitude
Source code at petsc4py/PETSc/DMPlex.pyx:2692
- Return type:
- metricGetMinimumMagnitude()#
Return the minimum tolerated metric magnitude.
Not collective.
See also
DMPlex.metricSetMinimumMagnitude
,DMPlex.metricGetMaximumMagnitude
,DMPlexMetricGetMinimumMagnitude
Source code at petsc4py/PETSc/DMPlex.pyx:2658
- Return type:
- metricGetNormalizationOrder()#
Return the order p for L-p normalization.
Not collective.
See also
DMPlex.metricSetNormalizationOrder
,DMPlex.metricGetTargetComplexity
,DMPlexMetricGetNormalizationOrder
Source code at petsc4py/PETSc/DMPlex.pyx:2794
- Return type:
- metricGetNumIterations()#
Return the number of parallel adaptation iterations.
Not collective.
Source code at petsc4py/PETSc/DMPlex.pyx:2624
- Return type:
- metricGetTargetComplexity()#
Return the target metric complexity.
Not collective.
See also
DMPlex.metricSetTargetComplexity
,DMPlex.metricGetNormalizationOrder
,DMPlexMetricGetTargetComplexity
Source code at petsc4py/PETSc/DMPlex.pyx:2760
- Return type:
- metricGetVerbosity()#
Return the verbosity of the mesh adaptation package.
Not collective.
- Returns:
verbosity – The verbosity, where -1 is silent and 10 is maximum.
- Return type:
- metricIntersection2(metric1, metric2, metricInt)#
Compute and return the intersection of two metrics.
Collective.
- Parameters:
- Return type:
- metricIntersection3(metric1, metric2, metric3, metricInt)#
Compute the intersection of three metrics.
Collective.
- Parameters:
- Return type:
- metricIsIsotropic()#
Return the flag indicating whether the metric is isotropic or not.
Not collective.
See also
DMPlex.metricSetIsotropic
,DMPlex.metricIsUniform
,DMPlex.metricRestrictAnisotropyFirst
,DMPlexMetricIsIsotropic
Source code at petsc4py/PETSc/DMPlex.pyx:2373
- Return type:
- metricIsUniform()#
Return the flag indicating whether the metric is uniform or not.
Not collective.
Source code at petsc4py/PETSc/DMPlex.pyx:2339
- Return type:
- metricNoInsertion()#
Return the flag indicating whether node insertion and deletion are turned off.
Not collective.
See also
DMPlex.metricSetNoInsertion
,DMPlex.metricNoSwapping
,DMPlex.metricNoMovement
,DMPlex.metricNoSurf
,DMPlexMetricNoInsertion
Source code at petsc4py/PETSc/DMPlex.pyx:2442
- Return type:
- metricNoMovement()#
Return the flag indicating whether node movement is turned off.
Not collective.
See also
DMPlex.metricSetNoMovement
,DMPlex.metricNoInsertion
,DMPlex.metricNoSwapping
,DMPlex.metricNoSurf
,DMPlexMetricNoMovement
Source code at petsc4py/PETSc/DMPlex.pyx:2514
- Return type:
- metricNoSurf()#
Return the flag indicating whether surface modification is turned off.
Not collective.
See also
DMPlex.metricSetNoSurf
,DMPlex.metricNoMovement
,DMPlex.metricNoInsertion
,DMPlex.metricNoSwapping
,DMPlexMetricNoSurf
Source code at petsc4py/PETSc/DMPlex.pyx:2550
- Return type:
- metricNoSwapping()#
Return the flag indicating whether facet swapping is turned off.
Not collective.
See also
DMPlex.metricSetNoSwapping
,DMPlex.metricNoInsertion
,DMPlex.metricNoMovement
,DMPlex.metricNoSurf
,DMPlexMetricNoSwapping
Source code at petsc4py/PETSc/DMPlex.pyx:2478
- Return type:
- metricNormalize(metric, ometric, determinant, restrictSizes=True, restrictAnisotropy=True)#
Apply L-p normalization to a metric.
Collective.
- Parameters:
metric (Vec) – The metric.
ometric (Vec) – The output metric.
determinant (Vec) – The output determinant.
restrictSizes (bool | None) – Flag indicating whether maximum/minimum magnitudes should be enforced.
restrictAnisotropy (bool | None) – Flag indicating whether maximum anisotropy should be enforced.
- Returns:
- Return type:
- metricRestrictAnisotropyFirst()#
Return
true
if anisotropy is restricted before normalization.Not collective.
See also
DMPlex.metricIsIsotropic
,DMPlex.metricSetRestrictAnisotropyFirst
,DMPlexMetricRestrictAnisotropyFirst
Source code at petsc4py/PETSc/DMPlex.pyx:2407
- Return type:
- metricSetFromOptions()#
Configure the object from the options database.
Collective.
See also
Source code at petsc4py/PETSc/DMPlex.pyx:2307
- Return type:
- metricSetGradationFactor(beta)#
Set the metric gradation factor.
Logically collective.
- metricSetHausdorffNumber(hausd)#
Set the metric Hausdorff number.
Logically collective.
- metricSetIsotropic(isotropic)#
Record whether the metric is isotropic or not.
Logically collective.
- Parameters:
isotropic (bool) – Flag indicating whether the metric is isotropic or not.
- Return type:
- metricSetMaximumAnisotropy(a_max)#
Set the maximum tolerated metric anisotropy.
Logically collective.
- metricSetMaximumMagnitude(h_max)#
Set the maximum tolerated metric magnitude.
Logically collective.
- metricSetMinimumMagnitude(h_min)#
Set the minimum tolerated metric magnitude.
Logically collective.
- metricSetNoInsertion(noInsert)#
Set the flag indicating whether node insertion should be turned off.
Logically collective.
- Parameters:
noInsert (bool) – Flag indicating whether node insertion and deletion should be turned off.
- Return type:
- metricSetNoMovement(noMove)#
Set the flag indicating whether node movement should be turned off.
Logically collective.
- Parameters:
noMove (bool) – Flag indicating whether node movement should be turned off.
- Return type:
- metricSetNoSurf(noSurf)#
Set the flag indicating whether surface modification should be turned off.
Logically collective.
- Parameters:
noSurf (bool) – Flag indicating whether surface modification should be turned off.
- Return type:
- metricSetNoSwapping(noSwap)#
Set the flag indicating whether facet swapping should be turned off.
Logically collective.
- Parameters:
noSwap (bool) – Flag indicating whether facet swapping should be turned off.
- Return type:
- metricSetNormalizationOrder(p)#
Set the order p for L-p normalization.
Logically collective.
- metricSetNumIterations(numIter)#
Set the number of parallel adaptation iterations.
Logically collective.
- metricSetRestrictAnisotropyFirst(restrictAnisotropyFirst)#
Record whether anisotropy is be restricted before normalization or after.
Logically collective.
- Parameters:
restrictAnisotropyFirst (bool) – Flag indicating if anisotropy is restricted before normalization or after.
- Return type:
- metricSetTargetComplexity(targetComplexity)#
Set the target metric complexity.
Logically collective.
- metricSetUniform(uniform)#
Record whether the metric is uniform or not.
Logically collective.
- Parameters:
uniform (bool) – Flag indicating whether the metric is uniform or not.
- Return type:
- metricSetVerbosity(verbosity)#
Set the verbosity of the mesh adaptation package.
Logically collective.
- Parameters:
verbosity (int) – The verbosity, where -1 is silent and 10 is maximum.
- Return type:
- orient()#
Give a consistent orientation to the input mesh.
Collective.
See also
Source code at petsc4py/PETSc/DMPlex.pyx:865
- Return type:
- permute(perm)#
Reorder the mesh according to the input permutation.
Collective.
- Parameters:
perm (IS) – The point permutation,
perm[old point number] = new point number
.- Returns:
pdm – The permuted
DMPlex
.- Return type:
See also
Redistribute shared points in order to achieve better balancing.
Collective.
- Parameters:
entityDepth (int | None) – Depth of the entity to balance (e.g., 0 -> balance vertices).
useInitialGuess (bool | None) – Whether to use the current distribution as initial guess.
parallel (bool | None) – Whether to use ParMETIS and do the partition in parallel or gather the graph onto a single process.
- Returns:
success – Whether the graph partitioning was successful or not. Unsuccessful simply means no change to the partitioning.
- Return type:
- reorderGetDefault()#
Return flag indicating whether the
DMPlex
should be reordered by default.Not collective.
See also
Source code at petsc4py/PETSc/DMPlex.pyx:2169
- Return type:
- reorderSetDefault(flag)#
Set flag indicating whether the DM should be reordered by default.
Logically collective.
- Parameters:
reorder – Flag for reordering.
flag (ReorderDefaultFlag)
- Return type:
- sectionLoad(viewer, sectiondm, sfxc)#
Load section into a
DM
.Collective.
- Parameters:
- Returns:
gsf (
SF
) – TheSF
that migrates any on-diskVec
data associated withsectionA
into a globalVec
associated with thesectiondm
’s global section (None
if not needed).lsf (
SF
) – TheSF
that migrates any on-diskVec
data associated withsectionA
into a localVec
associated with thesectiondm
’s local section (None
if not needed).
- Return type:
- sectionView(viewer, sectiondm)#
Save a section associated with a
DMPlex
.Collective.
- Parameters:
- Return type:
- setAdjacencyUseAnchors(useAnchors=True)#
Define adjacency in the mesh using the point-to-point constraints.
Logically collective.
- Parameters:
useAnchors (bool) – Flag to use the constraints. If
True
, then constrained points are omitted fromDMPlex.getAdjacency
, and their anchor points appear in their place.- Return type:
- setCellType(p, ctype)#
Set the polytope type of a given cell.
Not collective.
- Parameters:
p (int) – The cell.
ctype (PolytopeType) – The polytope type of the cell.
- Return type:
- setChart(pStart, pEnd)#
Set the interval for all mesh points [
pStart
,pEnd
).Not collective.
- Parameters:
- Return type:
See also
- setCone(p, cone, orientation=None)#
Set the points on the in-edges for this point in the DAG.
Not collective.
- Parameters:
- Return type:
- setConeOrientation(p, orientation)#
Set the orientations on the in-edges for this point in the DAG.
Not collective.
- Parameters:
p (int) – The point, which must lie in the chart set with
DMPlex.setChart
.
- Return type:
- setConeSize(p, size)#
Set the number of in-edges for this point in the DAG.
Not collective.
- Parameters:
p (int) – The point, which must lie in the chart set with
DMPlex.setChart
.size (int) – The cone size for point
p
.
- Return type:
See also
DM
,DMPlex
,DMPlex.create
,DMPlex.getConeSize
,DMPlex.setChart
,DMPlexSetConeSize
- setMatClosure(sec, gsec, mat, point, values, addv=None)#
Set an array of the values on the closure of
point
.Not collective.
- Parameters:
sec (Section) – The section describing the layout in
mat
, orNone
to use the default section.gsec (Section) – The section describing the layout in
mat
, orNone
to use the default global section.mat (Mat) – The matrix.
values (Sequence[Scalar]) – The array of values.
mode – The insertion mode.
addv (InsertModeSpec | None)
- Return type:
See also
- setPartitioner(part)#
Set the mesh partitioner.
Logically collective.
- Parameters:
part (Partitioner) – The partitioner.
- Return type:
- setRefinementLimit(refinementLimit)#
Set the maximum cell volume for refinement.
Logically collective.
- Parameters:
refinementLimit (float) – The maximum cell volume in the refined mesh.
- Return type:
- setRefinementUniform(refinementUniform=True)#
Set the flag for uniform refinement.
Logically collective.
- setSupport(p, supp)#
Set the points on the out-edges for this point in the DAG.
Not collective.
- Parameters:
p (int) – The point, which must lie in the chart set with
DMPlex.setChart
.supp (Sequence[int]) – An array of points which are on the out-edges for point
p
.
- Return type:
- setSupportSize(p, size)#
Set the number of out-edges for this point in the DAG.
Not collective.
- Parameters:
p (int) – The point, which must lie in the chart set with
DMPlex.setChart
.size (int) – The support size for point
p
.
- Return type:
- setTetGenOptions(opts)#
Set the options used for the Tetgen mesh generator.
Not collective.
- setTriangleOptions(opts)#
Set the options used for the Triangle mesh generator.
Not collective.
- setVecClosure(sec, vec, point, values, addv=None)#
Set an array of the values on the closure of
point
.Not collective.
- Parameters:
- Return type:
See also
- stratify()#
Calculate the strata of DAG.
Collective.
See also
DM
,DMPlex
,DMPlex.create
,DMPlex.symmetrize
,DMPlexStratify
Source code at petsc4py/PETSc/DMPlex.pyx:853
- Return type:
- symmetrize()#
Create support (out-edge) information from cone (in-edge) information.
Not collective.
See also
DM
,DMPlex
,DMPlex.create
,DMPlex.setChart
,DMPlex.setConeSize
,DMPlex.setCone
,DMPlexSymmetrize
Source code at petsc4py/PETSc/DMPlex.pyx:840
- Return type:
- topologyLoad(viewer)#
Load a topology into this
DMPlex
object.Collective.
- Parameters:
- Returns:
sfxc – The
SF
that pushes points in[0, N)
to the associated points in the loadedDMPlex
, whereN
is the global number of points.- Return type:
See also
DM
,DMPlex
,DM.load
,DMPlex.coordinatesLoad
,DMPlex.labelsLoad
,DM.view
,SF
,Viewer
,DMPlexTopologyLoad
- uninterpolate()#
Convert to a mesh with only cells and vertices.
Collective.
Source code at petsc4py/PETSc/DMPlex.pyx:1730
- Return type:
- vecGetClosure(sec, vec, p)#
Return an array of values on the closure of
p
.Not collective.
- Parameters:
- Return type:
See also