petsc4py.PETSc.PC#
- class petsc4py.PETSc.PC#
Bases:
Object
Preconditioners.
PC
is described in thePETSc manual
. Calling thePC
with a vector as an argument willapply
the preconditioner as shown in the example below.Examples
>>> from petsc4py import PETSc >>> v = PETSc.Vec().createWithArray([1, 2]) >>> m = PETSc.Mat().createDense(2, array=[[1, 0], [0, 1]]) >>> pc = PETSc.PC().create() >>> pc.setOperators(m) >>> u = pc(v) # u is created internally >>> pc.apply(v, u) # u can also be passed as second argument
See also
Enumerations
The ASM subtype.
The composite type.
The deflation space subtype.
The reason the preconditioner has failed.
The field split Schur factorization type.
The field split Schur subtype.
The GAMG subtype.
The GASM subtype.
The HPDDM coarse correction type.
The MG cycle type.
The MG subtype.
The patch construction type.
The manner in which the preconditioner is applied.
The preconditioner method.
Methods Summary
addCompositePCType
(pc_type)appendOptionsPrefix
(prefix)Append to the prefix used for all the
PC
options.apply
(x, y)Apply the
PC
to a vector.applySymmetricLeft
(x, y)Apply the left part of a symmetric
PC
to a vector.applySymmetricRight
(x, y)Apply the right part of a symmetric
PC
to a vector.applyTranspose
(x, y)Apply the transpose of the
PC
to a vector.create
([comm])Create an empty
PC
.createPython
([context, comm])Create a preconditioner of Python type.
destroy
()Return the local
KSP
object for all blocks on this process.Return a component of the composite
PC
.getDM
()Return the coarse problem
KSP
.Return the additional preconditioner.
Return the factored matrix.
Return the solver package used to perform the factorization.
Return the reason the
PC
terminated.Return the
KSP
for the Schur complement based splits.getFieldSplitSubIS
(splitname)Return the
IS
associated with a given name.Return the
KSP
for all splits.Return the coarse correction type.
Compute the grid and operator complexities.
Return true if the
KSP
in SLEPcST
and the subdomain solver is shared.Return the
Type.HYPRE
type.getKSP
()Return the
KSP
used on the coarse grid.getMGInterpolation
(level)Return the interpolation operator for the given level.
Return the number of
MG
levels.getMGRScale
(level)Return the pointwise scaling for the restriction operator on the given level.
getMGRestriction
(level)Return the restriction operator for the given level.
getMGSmoother
(level)Return the
KSP
to be used as a smoother.getMGSmootherDown
(level)Return the
KSP
to be used as a smoother before coarse grid correction.getMGSmootherUp
(level)Return the
KSP
to be used as a smoother after coarse grid correction.Return the form of multigrid.
Return the matrices associated with a linear system.
Return the prefix used for all the
PC
options.Return the local
KSP
object for all blocks on this process.Return the instance of the class implementing the required Python methods.
Return the fully qualified Python name of the class used by the preconditioner.
getType
()Return the preconditioner type.
Return the flag to indicate if
PC
is applied toA
orP
.matApply
(x, y)Apply the
PC
to many vectors stored asMat.Type.DENSE
.reset
()Reset the
PC
, removing any allocated vectors and matrices.setASMLocalSubdomains
(nsd[, is_sub, is_local])Set the local subdomains.
setASMOverlap
(overlap)Set the overlap between a pair of subdomains.
setASMSortIndices
(dosort)Set to sort subdomain indices.
setASMTotalSubdomains
(nsd[, is_sub, is_local])Set the subdomains for all processes.
setASMType
(asmtype)Set the type of restriction and interpolation.
setBDDCChangeOfBasisMat
(T[, interior])Set a user defined change of basis for degrees of freedom.
setBDDCCoarseningRatio
(cratio)Set the coarsening ratio used in the multilevel version.
Set the
IS
defining Dirichlet boundaries for the global problem.Set the
IS
defining Dirichlet boundaries in local ordering.setBDDCDiscreteGradient
(G[, order, field, ...])Set the discrete gradient.
setBDDCDivergenceMat
(div[, trans, l2l])Set the linear operator representing ∫ div(u)•p dx.
setBDDCDofsSplitting
(isfields)Set the index set(s) defining fields of the global matrix.
setBDDCDofsSplittingLocal
(isfields)Set the index set(s) defining fields of the local subdomain matrix.
setBDDCLevels
(levels)Set the maximum number of additional levels allowed.
Provide a custom connectivity graph for local dofs.
setBDDCNeumannBoundaries
(bndr)Set the
IS
defining Neumann boundaries for the global problem.Set the
IS
defining Neumann boundaries in local ordering.setBDDCPrimalVerticesIS
(primv)Set additional user defined primal vertices.
setBDDCPrimalVerticesLocalIS
(primv)Set additional user defined primal vertices.
setCompositeType
(ctype)Set the type of composite preconditioner.
setCoordinates
(coordinates)Set the coordinates for the nodes on the local process.
setDM
(dm)Set the
DM
that may be used by some preconditioners.Set the coarse problem matrix.
Set the coarse problem correction factor.
setDeflationInitOnly
(flg)Set to only perform the initialization.
setDeflationLevels
(levels)Set the maximum level of deflation nesting.
Set the projection null space matrix.
Set the reduction factor for the preconditioner.
setDeflationSpace
(W, transpose)Set the deflation space matrix or its (Hermitian) transpose.
setDeflationSpaceToCompute
(space_type, size)Set the deflation space type.
setFactorLevels
(levels)Set the number of levels of fill.
setFactorOrdering
([ord_type, nzdiag, reuse])Set options for the matrix factorization reordering.
setFactorPivot
([zeropivot, inblocks])Set options for matrix factorization pivoting.
Set up the factorization solver.
setFactorShift
([shift_type, amount])Set options for shifting diagonal entries of a matrix.
setFactorSolverType
(solver)Set the solver package used to perform the factorization.
setFailedReason
(reason)Set the reason the
PC
terminated.setFieldSplitFields
(bsize, *fields)Sets the elements for the field split.
setFieldSplitIS
(*fields)Set the elements for the field split by
IS
.setFieldSplitSchurFactType
(ctype)Set the type of approximate block factorization.
setFieldSplitSchurPreType
(ptype[, pre])Set from what operator the
PC
is constructed.setFieldSplitType
(ctype)Set the type of composition of a field split preconditioner.
Set various
PC
parameters from user options.setGAMGLevels
(levels)Set the maximum number of levels.
setGAMGSmooths
(smooths)Set the number of smoothing steps used on all levels.
setGAMGType
(gamgtype)Set the type of algorithm.
setGASMOverlap
(overlap)Set the overlap between a pair of subdomains.
setGASMType
(gasmtype)Set the type of restriction and interpolation.
setHPDDMAuxiliaryMat
(uis, uaux)Set the auxiliary matrix used by the preconditioner.
setHPDDMCoarseCorrectionType
(correction_type)Set the coarse correction type.
setHPDDMDeflationMat
(uis, U)Set the deflation space used to assemble a coarse operator.
Set to indicate that the
Mat
passed to thePC
is the local Neumann matrix.Set the right-hand side matrix of the preconditioner.
setHYPREAMSSetInteriorNodes
(interior)Set the list of interior nodes to a zero conductivity region.
setHYPREDiscreteCurl
(mat)Set the discrete curl matrix.
Set the discrete gradient matrix.
Set the vector Poisson matrix.
setHYPRESetBetaPoissonMatrix
([mat])Set the Posson matrix.
setHYPRESetEdgeConstantVectors
(ozz, zoz[, zzo])Set the representation of the constant vector fields in the edge element basis.
setHYPRESetInterpolations
(dim[, RT_Pi_Full, ...])Set the interpolation matrices.
setHYPREType
(hypretype)Set the
Type.HYPRE
type.setMGCycleType
(cycle_type)Set the type of cycles.
setMGCycleTypeOnLevel
(level, cycle_type)Set the type of cycle on the given level.
setMGInterpolation
(level, mat)Set the interpolation operator for the given level.
setMGLevels
(levels)Set the number of
MG
levels.setMGR
(level, r)Set the vector where the residual is stored.
setMGRScale
(level, rscale)Set the pointwise scaling for the restriction operator on the given level.
setMGRestriction
(level, mat)Set the restriction operator for the given level.
setMGRhs
(level, rhs)Set the vector where the right-hand side is stored.
setMGType
(mgtype)Set the form of multigrid.
setMGX
(level, x)Set the vector where the solution is stored.
setOperators
([A, P])Set the matrices associated with the linear system.
setOptionsPrefix
(prefix)Set the prefix used for all the
PC
options.Set the cell numbering.
setPatchComputeFunction
(function[, args, kargs])Set compute operator callbacks.
setPatchComputeFunctionInteriorFacets
(function)Set compute operator callbacks.
setPatchComputeOperator
(operator[, args, kargs])Set compute operator callbacks.
setPatchComputeOperatorInteriorFacets
(operator)Set compute operator callbacks.
setPatchConstructType
(typ[, operator, args, ...])Set compute operator callbacks.
setPatchDiscretisationInfo
(dms, bs, ...)Set discretisation info.
setPythonContext
(context)Set the instance of the class implementing the required Python methods.
setPythonType
(py_type)Set the fully qualified Python name of the class to be used.
setReusePreconditioner
(flag)Set to indicate the preconditioner is to be reused.
Set the block size of the preconditioner.
setSPAICacheSize
(size)Set the cache size.
setSPAIEpsilon
(val)Set the tolerance for the preconditioner.
setSPAIMax
(maxval)Set the size of working buffers in the preconditioner.
setSPAIMaxNew
(maxval)Set the maximum number of new non-zero candidates per step.
setSPAINBSteps
(nbsteps)Set the maximum number of improvement steps per row.
setSPAISp
(sym)Set to specify a symmetric sparsity pattern.
setSPAIVerbose
(level)Set the verbosity level.
setType
(pc_type)Set the preconditioner type.
setUp
()Set up the internal data structures for the
PC
.Set up the
PC
for each block.setUseAmat
(flag)Set to indicate to apply
PC
toA
and notP
.view
([viewer])View the
PC
object.Methods Documentation
- appendOptionsPrefix(prefix)#
Append to the prefix used for all the
PC
options.Logically collective.
- applySymmetricLeft(x, y)#
Apply the left part of a symmetric
PC
to a vector.Collective.
- Parameters:
- Return type:
See also
- applySymmetricRight(x, y)#
Apply the right part of a symmetric
PC
to a vector.Collective.
- Parameters:
- Return type:
See also
- applyTranspose(x, y)#
Apply the transpose of the
PC
to a vector.Collective.
For complex numbers this applies the non-Hermitian transpose.
- Parameters:
- Return type:
See also
- create(comm=None)#
Create an empty
PC
.Collective.
The default preconditioner for sparse matrices is
ILU
orICC
with 0 fill on one process and block Jacobi (BJACOBI
) withILU
orICC
in parallel. For dense matrices it is alwaysNone
.- Parameters:
comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.- Return type:
- createPython(context=None, comm=None)#
Create a preconditioner of Python type.
Collective.
- Parameters:
context (Any) – An instance of the Python class implementing the required methods.
comm (Comm | None) – MPI communicator, defaults to
Sys.getDefaultComm
.
- Return type:
- destroy()#
Destroy the
PC
that was created withcreate
.Collective.
See also
Source code at petsc4py/PETSc/PC.pyx:250
- Return type:
- getASMSubKSP()#
Return the local
KSP
object for all blocks on this process.Not collective.
See also
- getDM()#
Return the
DM
associated with thePC
.Not collective.
See also
Source code at petsc4py/PETSc/PC.pyx:681
- Return type:
- getDeflationCoarseKSP()#
Return the coarse problem
KSP
.Not collective.
See also
Source code at petsc4py/PETSc/PC.pyx:2910
- Return type:
- getDeflationPC()#
Return the additional preconditioner.
Not collective.
See also
Source code at petsc4py/PETSc/PC.pyx:2925
- Return type:
- getFactorMatrix()#
Return the factored matrix.
Not collective.
See also
Source code at petsc4py/PETSc/PC.pyx:1434
- Return type:
- getFactorSolverType()#
Return the solver package used to perform the factorization.
Not collective.
See also
Source code at petsc4py/PETSc/PC.pyx:1290
- Return type:
- getFailedReason()#
Return the reason the
PC
terminated.Not collective.
After a call to KSPCheckDot() or KSPCheckNorm() inside a KSPSolve(), or after a call to PCReduceFailedReason() this is the maximum reason over all ranks in the PC communicator and hence logically collective. Otherwise it is the local value.
See also
Source code at petsc4py/PETSc/PC.pyx:520
- Return type:
- getFieldSplitSchurGetSubKSP()#
Return the
KSP
for the Schur complement based splits.Not collective.
- getHPDDMCoarseCorrectionType()#
Return the coarse correction type.
Not collective.
See also
Source code at petsc4py/PETSc/PC.pyx:2562
- Return type:
- getHPDDMComplexities()#
Compute the grid and operator complexities.
Collective.
See also
Return true if the
KSP
in SLEPcST
and the subdomain solver is shared.Not collective.
See also
Source code at petsc4py/PETSc/PC.pyx:2576
- Return type:
- getHYPREType()#
Return the
Type.HYPRE
type.Not collective.
See also
Source code at petsc4py/PETSc/PC.pyx:1076
- Return type:
- getKSP()#
Return the
KSP
if thePC
isType.KSP
.Not collective.
See also
Source code at petsc4py/PETSc/PC.pyx:1686
- Return type:
- getMGCoarseSolve()#
Return the
KSP
used on the coarse grid.Not collective.
See also
Source code at petsc4py/PETSc/PC.pyx:1762
- Return type:
- getMGInterpolation(level)#
Return the interpolation operator for the given level.
Logically collective.
- Parameters:
level (int) – The level where interpolation is defined from
level-1
tolevel
.- Return type:
See also
- getMGLevels()#
Return the number of
MG
levels.Not collective.
See also
Source code at petsc4py/PETSc/PC.pyx:1730
- Return type:
- getMGRScale(level)#
Return the pointwise scaling for the restriction operator on the given level.
Logically collective.
- Parameters:
level (int) – The level where restriction is defined from
level
tolevel-1
.- Return type:
See also
- getMGRestriction(level)#
Return the restriction operator for the given level.
Logically collective.
- Parameters:
level (int) – The level where restriction is defined from
level
tolevel-1
.- Return type:
See also
- getMGSmootherDown(level)#
Return the
KSP
to be used as a smoother before coarse grid correction.Not collective.
See also
- getMGSmootherUp(level)#
Return the
KSP
to be used as a smoother after coarse grid correction.Not collective.
See also
- getMGType()#
Return the form of multigrid.
Logically collective.
See also
Source code at petsc4py/PETSc/PC.pyx:1703
- Return type:
- getOperators()#
Return the matrices associated with a linear system.
Not collective.
See also
- getOptionsPrefix()#
Return the prefix used for all the
PC
options.Not collective.
See also
Source code at petsc4py/PETSc/PC.pyx:341
- Return type:
- getPythonContext()#
Return the instance of the class implementing the required Python methods.
Not collective.
Source code at petsc4py/PETSc/PC.pyx:779
- Return type:
- getPythonType()#
Return the fully qualified Python name of the class used by the preconditioner.
Not collective.
Source code at petsc4py/PETSc/PC.pyx:808
- Return type:
- getType()#
Return the preconditioner type.
Not collective.
Source code at petsc4py/PETSc/PC.pyx:308
- Return type:
- getUseAmat()#
Return the flag to indicate if
PC
is applied toA
orP
.Logically collective.
- Returns:
flag – True if
A
is used and False ifP
.- Return type:
See also
- matApply(x, y)#
Apply the
PC
to many vectors stored asMat.Type.DENSE
.Collective.
- Parameters:
- Return type:
See also
- reset()#
Reset the
PC
, removing any allocated vectors and matrices.Collective.
See also
Source code at petsc4py/PETSc/PC.pyx:553
- Return type:
- setASMLocalSubdomains(nsd, is_sub=None, is_local=None)#
Set the local subdomains.
Collective.
- Parameters:
- Return type:
See also
- setASMOverlap(overlap)#
Set the overlap between a pair of subdomains.
Logically collective.
See also
- setASMSortIndices(dosort)#
Set to sort subdomain indices.
Logically collective.
See also
- setASMTotalSubdomains(nsd, is_sub=None, is_local=None)#
Set the subdomains for all processes.
Collective.
- Parameters:
- Return type:
See also
- setASMType(asmtype)#
Set the type of restriction and interpolation.
Logically collective.
See also
- setBDDCChangeOfBasisMat(T, interior=False)#
Set a user defined change of basis for degrees of freedom.
Collective.
- Parameters:
- Return type:
See also
- setBDDCCoarseningRatio(cratio)#
Set the coarsening ratio used in the multilevel version.
Logically collective.
See also
- setBDDCDirichletBoundaries(bndr)#
Set the
IS
defining Dirichlet boundaries for the global problem.Collective.
See also
- setBDDCDirichletBoundariesLocal(bndr)#
Set the
IS
defining Dirichlet boundaries in local ordering.Collective.
- setBDDCDiscreteGradient(G, order=1, field=1, gord=True, conforming=True)#
Set the discrete gradient.
Collective.
- Parameters:
G (Mat) – The discrete gradient matrix in
Mat.Type.AIJ
format.order (int) – The order of the Nedelec space.
field (int) – The field number of the Nedelec degrees of freedom. This is not used if no fields have been specified.
gord (bool) – Enable to use global ordering in the rows of
G
.conforming (bool) – Enable if the mesh is conforming.
- Return type:
See also
- setBDDCDivergenceMat(div, trans=False, l2l=None)#
Set the linear operator representing ∫ div(u)•p dx.
Collective.
- Parameters:
- Return type:
See also
- setBDDCDofsSplitting(isfields)#
Set the index set(s) defining fields of the global matrix.
Collective.
- Parameters:
isfields (IS | Sequence[IS]) – The sequence of
IS
describing the fields in global ordering.- Return type:
See also
- setBDDCDofsSplittingLocal(isfields)#
Set the index set(s) defining fields of the local subdomain matrix.
Collective.
Not all nodes need to be listed. Unlisted nodes will belong to the complement field.
- Parameters:
isfields (IS | Sequence[IS]) – The sequence of
IS
describing the fields in local ordering.- Return type:
See also
- setBDDCLevels(levels)#
Set the maximum number of additional levels allowed.
Logically collective.
See also
- setBDDCLocalAdjacency(csr)#
Provide a custom connectivity graph for local dofs.
Not collective.
- Parameters:
csr (CSRIndicesSpec) – Compressed sparse row layout information.
- Return type:
See also
- setBDDCNeumannBoundaries(bndr)#
Set the
IS
defining Neumann boundaries for the global problem.Collective.
See also
- setBDDCNeumannBoundariesLocal(bndr)#
Set the
IS
defining Neumann boundaries in local ordering.Collective.
- setBDDCPrimalVerticesIS(primv)#
Set additional user defined primal vertices.
Collective.
See also
- setBDDCPrimalVerticesLocalIS(primv)#
Set additional user defined primal vertices.
Collective.
See also
- setCompositeType(ctype)#
Set the type of composite preconditioner.
Logically collective.
- Parameters:
ctype (CompositeType) – The type of composition.
- Return type:
See also
- setCoordinates(coordinates)#
Set the coordinates for the nodes on the local process.
Collective.
- Parameters:
coordinates (Sequence[Sequence[float]]) – The two dimensional coordinate array.
- Return type:
See also
- setDeflationCoarseMat(mat)#
Set the coarse problem matrix.
Collective.
See also
- setDeflationCorrectionFactor(fact)#
Set the coarse problem correction factor.
Logically collective.
See also
- setDeflationInitOnly(flg)#
Set to only perform the initialization.
Logically collective.
Sets initial guess to the solution on the deflation space but does not apply the deflation preconditioner. The additional preconditioner is still applied.
See also
- setDeflationLevels(levels)#
Set the maximum level of deflation nesting.
Logically collective.
See also
- setDeflationProjectionNullSpaceMat(mat)#
Set the projection null space matrix.
Collective.
See also
- setDeflationReductionFactor(red)#
Set the reduction factor for the preconditioner.
Logically collective.
See also
- setDeflationSpace(W, transpose)#
Set the deflation space matrix or its (Hermitian) transpose.
Logically collective.
- Parameters:
- Return type:
See also
- setDeflationSpaceToCompute(space_type, size)#
Set the deflation space type.
Logically collective.
- Parameters:
space_type (DeflationSpaceType) – The deflation space type.
size (int) – The size of the space to compute
- Return type:
See also
- setFactorLevels(levels)#
Set the number of levels of fill.
Logically collective.
See also
- setFactorOrdering(ord_type=None, nzdiag=None, reuse=None)#
Set options for the matrix factorization reordering.
Logically collective.
- Parameters:
- Return type:
- setFactorPivot(zeropivot=None, inblocks=None)#
Set options for matrix factorization pivoting.
Logically collective.
- Parameters:
- Return type:
See also
- setFactorSetUpSolverType()#
Set up the factorization solver.
Collective.
This can be called after
KSP.setOperators
orPC.setOperators
, causesMatGetFactor
to be called so then one may set the options for that particular factorization object.Source code at petsc4py/PETSc/PC.pyx:1304
- Return type:
- setFactorShift(shift_type=None, amount=None)#
Set options for shifting diagonal entries of a matrix.
Logically collective.
- Parameters:
- Return type:
See also
- setFactorSolverType(solver)#
Set the solver package used to perform the factorization.
Logically collective.
- Parameters:
solver (SolverType | str) – The solver package used to factorize.
- Return type:
See also
- setFailedReason(reason)#
Set the reason the
PC
terminated.Logically collective.
- Parameters:
reason (FailedReason | str) – the reason the
PC
terminated- Return type:
See also
- setFieldSplitFields(bsize, *fields)#
Sets the elements for the field split.
Collective.
- Parameters:
- Return type:
- setFieldSplitIS(*fields)#
Set the elements for the field split by
IS
.Logically collective.
Solve options for this split will be available under the prefix
-fieldsplit_SPLITNAME_*
.- Parameters:
fields (Tuple[str, IS]) – A sequence of tuples containing the split name and the
IS
that defines the elements in the split.- Return type:
See also
- setFieldSplitSchurFactType(ctype)#
Set the type of approximate block factorization.
Collective.
- Parameters:
ctype (FieldSplitSchurFactType) – The type indicating which blocks to retain.
- Return type:
See also
- setFieldSplitSchurPreType(ptype, pre=None)#
Set from what operator the
PC
is constructed.Collective.
- Parameters:
ptype (FieldSplitSchurPreType) – The type of matrix to use for preconditioning the Schur complement.
pre (Mat | None) – The optional matrix to use for preconditioning.
- Return type:
See also
- setFieldSplitType(ctype)#
Set the type of composition of a field split preconditioner.
Collective.
- Parameters:
ctype (CompositeType) – The type of composition.
- Return type:
See also
- setFromOptions()#
Set various
PC
parameters from user options.Collective.
See also
Source code at petsc4py/PETSc/PC.pyx:374
- Return type:
- setGAMGLevels(levels)#
Set the maximum number of levels.
Not collective.
See also
- setGAMGSmooths(smooths)#
Set the number of smoothing steps used on all levels.
Logically collective.
See also
- setGAMGType(gamgtype)#
Set the type of algorithm.
Collective.
See also
- setGASMOverlap(overlap)#
Set the overlap between a pair of subdomains.
Logically collective.
See also
- setGASMType(gasmtype)#
Set the type of restriction and interpolation.
Logically collective.
See also
- setHPDDMAuxiliaryMat(uis, uaux)#
Set the auxiliary matrix used by the preconditioner.
Logically collective.
- Parameters:
- Return type:
See also
- setHPDDMCoarseCorrectionType(correction_type)#
Set the coarse correction type.
Collective.
- Parameters:
correction_type (HPDDMCoarseCorrectionType) – The type of coarse correction to apply.
- Return type:
See also
- setHPDDMDeflationMat(uis, U)#
Set the deflation space used to assemble a coarse operator.
Logically collective.
- Parameters:
U (Mat) – The deflation sequential matrix of type
Mat.Type.DENSE
.
- Return type:
See also
- setHPDDMHasNeumannMat(has)#
Set to indicate that the
Mat
passed to thePC
is the local Neumann matrix.Logically collective.
- Parameters:
has (bool) – Enable to indicate the matrix is the local Neumann matrix.
- Return type:
See also
- setHPDDMRHSMat(B)#
Set the right-hand side matrix of the preconditioner.
Logically collective.
See also
- setHYPREAMSSetInteriorNodes(interior)#
Set the list of interior nodes to a zero conductivity region.
Collective.
- Parameters:
interior (Vec) – A vector where a value of 1.0 indicates an interior node.
- Return type:
See also
- setHYPREDiscreteCurl(mat)#
Set the discrete curl matrix.
Collective.
See also
- setHYPREDiscreteGradient(mat)#
Set the discrete gradient matrix.
Collective.
See also
- setHYPRESetAlphaPoissonMatrix(mat)#
Set the vector Poisson matrix.
Collective.
See also
- setHYPRESetBetaPoissonMatrix(mat=None)#
Set the Posson matrix.
Collective.
See also
- setHYPRESetEdgeConstantVectors(ozz, zoz, zzo=None)#
Set the representation of the constant vector fields in the edge element basis.
Collective.
- Parameters:
- Return type:
See also
- setHYPRESetInterpolations(dim, RT_Pi_Full=None, RT_Pi=None, ND_Pi_Full=None, ND_Pi=None)#
Set the interpolation matrices.
Collective.
- Parameters:
dim (int) – The dimension of the problem.
RT_Pi_Full (Mat | None) – The Raviart-Thomas interpolation matrix or
None
to omit.RT_Pi – The xyz components of the Raviart-Thomas interpolation matrix, or
None
to omit.ND_Pi_Full (Mat | None) – The Nedelec interpolation matrix or
None
to omit.ND_Pi – The xyz components of the Nedelec interpolation matrix, or
None
to omit.
- Return type:
See also
- setHYPREType(hypretype)#
Set the
Type.HYPRE
type.Collective.
- Parameters:
hypretype (str) – The name of the type, one of
"euclid"
,"pilut"
,"parasails"
,"boomeramg"
,"ams"
,"ads"
- Return type:
See also
- setMGCycleType(cycle_type)#
Set the type of cycles.
Logically collective.
- Parameters:
cycle_type (MGCycleType) – The type of multigrid cycles to use.
- Return type:
See also
- setMGCycleTypeOnLevel(level, cycle_type)#
Set the type of cycle on the given level.
Logically collective.
- Parameters:
level (int) – The level on which to set the cycle type.
cycle_type (MGCycleType) – The type of multigrid cycles to use.
- Return type:
See also
- setMGInterpolation(level, mat)#
Set the interpolation operator for the given level.
Logically collective.
- Parameters:
level – The level where interpolation is defined from
level-1
tolevel
.mat (Mat) – The interpolation operator
- Return type:
See also
- setMGR(level, r)#
Set the vector where the residual is stored.
Logically collective.
If not provided, one will be set internally. Will be cleaned up in
destroy
.- Parameters:
- Return type:
See also
- setMGRScale(level, rscale)#
Set the pointwise scaling for the restriction operator on the given level.
Logically collective.
- Parameters:
- Return type:
See also
- setMGRestriction(level, mat)#
Set the restriction operator for the given level.
Logically collective.
- Parameters:
- Return type:
See also
- setMGRhs(level, rhs)#
Set the vector where the right-hand side is stored.
Logically collective.
If not provided, one will be set internally. Will be cleaned up in
destroy
.- Parameters:
- Return type:
See also
- setMGType(mgtype)#
Set the form of multigrid.
Logically collective.
See also
- setMGX(level, x)#
Set the vector where the solution is stored.
Logically collective.
If not provided, one will be set internally. Will be cleaned up in
destroy
.- Parameters:
- Return type:
See also
- setOperators(A=None, P=None)#
Set the matrices associated with the linear system.
Logically collective.
Passing
None
forA
orP
removes the matrix that is currently used. PETSc does not reset the matrix entries of eitherA
orP
to zero after a linear solve; the user is completely responsible for matrix assembly. SeeMat.zeroEntries
to zero all elements of a matrix.- Parameters:
- Return type:
See also
- setOptionsPrefix(prefix)#
Set the prefix used for all the
PC
options.Logically collective.
See also
- setPatchCellNumbering(sec)#
Set the cell numbering.
- setPatchComputeFunction(function, args=None, kargs=None)#
Set compute operator callbacks.
Source code at petsc4py/PETSc/PC.pyx:2444
- Return type:
- setPatchComputeFunctionInteriorFacets(function, args=None, kargs=None)#
Set compute operator callbacks.
Source code at petsc4py/PETSc/PC.pyx:2452
- Return type:
- setPatchComputeOperator(operator, args=None, kargs=None)#
Set compute operator callbacks.
Source code at petsc4py/PETSc/PC.pyx:2428
- Return type:
- setPatchComputeOperatorInteriorFacets(operator, args=None, kargs=None)#
Set compute operator callbacks.
Source code at petsc4py/PETSc/PC.pyx:2436
- Return type:
- setPatchConstructType(typ, operator=None, args=None, kargs=None)#
Set compute operator callbacks.
Source code at petsc4py/PETSc/PC.pyx:2460
- Return type:
- setPatchDiscretisationInfo(dms, bs, cellNodeMaps, subspaceOffsets, ghostBcNodes, globalBcNodes)#
Set discretisation info.
Source code at petsc4py/PETSc/PC.pyx:2387
- Return type:
- setPythonContext(context)#
Set the instance of the class implementing the required Python methods.
Not collective.
- setPythonType(py_type)#
Set the fully qualified Python name of the class to be used.
Collective.
- setReusePreconditioner(flag)#
Set to indicate the preconditioner is to be reused.
Logically collective.
Normally if the
A
matrix inside aPC
changes, thePC
automatically updates itself using information from the changed matrix. Enable this option prevents this.- Parameters:
flag (bool) – Set to
True
to use the reuse the current preconditioner andFalse
to recompute on changes to the matrix.- Return type:
See also
- setSPAIBlockSize(n)#
Set the block size of the preconditioner.
Logically collective.
See also
- setSPAICacheSize(size)#
Set the cache size.
Logically collective.
See also
- setSPAIEpsilon(val)#
Set the tolerance for the preconditioner.
Logically collective.
See also
- setSPAIMax(maxval)#
Set the size of working buffers in the preconditioner.
Logically collective.
- Parameters:
maxval (int) – Number of entries in the work arrays to be allocated, defaults to
5000
.- Return type:
See also
- setSPAIMaxNew(maxval)#
Set the maximum number of new non-zero candidates per step.
Logically collective.
See also
- setSPAINBSteps(nbsteps)#
Set the maximum number of improvement steps per row.
Logically collective.
See also
- setSPAISp(sym)#
Set to specify a symmetric sparsity pattern.
Logically collective.
See also
- setSPAIVerbose(level)#
Set the verbosity level.
Logically collective.
See also
- setType(pc_type)#
Set the preconditioner type.
Collective.
See also
- setUp()#
Set up the internal data structures for the
PC
.Collective.
See also
Source code at petsc4py/PETSc/PC.pyx:541
- Return type:
- setUpOnBlocks()#
Set up the
PC
for each block.Collective.
For nested preconditioners such as
BJACOBI
,setUp
is not called on each sub-KSP
whensetUp
is called on the outerPC
. This routine ensures it is called.See also
Source code at petsc4py/PETSc/PC.pyx:565
- Return type:
- setUseAmat(flag)#
Set to indicate to apply
PC
toA
and notP
.Logically collective.
Sets a flag to indicate that when the preconditioner needs to apply (part of) the operator during the preconditioning process, it applies to
A
provided toTS.setRHSJacobian
,TS.setIJacobian
,SNES.setJacobian
,KSP.setOperators
orPC.setOperators
not theP
.See also