petsc4py.PETSc.PC#

class petsc4py.PETSc.PC#

Bases: Object

Preconditioners.

PC is described in the PETSc manual. Calling the PC with a vector as an argument will apply 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) # Vec u is created internally, can also be passed as second argument

See also

PC

Enumerations

ASMType

The ASM subtype.

CompositeType

The composite type.

DeflationSpaceType

The deflation space subtype.

FailedReason

The reason the preconditioner has failed.

FieldSplitSchurFactType

The field split Schur factorization type.

FieldSplitSchurPreType

The field split Schur subtype.

GAMGType

The GAMG subtype.

GASMType

The GASM subtype.

HPDDMCoarseCorrectionType

The HPDDM coarse correction type.

MGCycleType

The MG cycle type.

MGType

The MG subtype.

PatchConstructType

The patch construction type.

Side

The manner in which the preconditioner is applied.

Type

The preconditioner method.

Methods Summary

addCompositePCType(pc_type)

Add a PC of the given type to the composite PC.

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()

Destroy the PC that was created with create.

getASMSubKSP()

Return the local KSP object for all blocks on this process.

getCompositePC(n)

Return a component of the composite PC.

getDM()

Return the DM associated with the PC.

getDeflationCoarseKSP()

Return the coarse problem KSP.

getDeflationPC()

Return the additional preconditioner.

getFactorMatrix()

Return the factored matrix.

getFactorSolverType()

Return the solver package used to perform the factorization.

getFailedReason()

Return the reason the PC terminated.

getFailedReasonRank()

Return the reason the PC terminated on this rank.

getFieldSplitSchurGetSubKSP()

Return the KSP for the Schur complement based splits.

getFieldSplitSubKSP()

Return the KSP for all splits.

getHPDDMCoarseCorrectionType()

Return the coarse correction type.

getHPDDMSTShareSubKSP()

Return true if the KSP in SLEPc ST and the subdomain solver is shared.

getHYPREType()

Return the Type.HYPRE type.

getKSP()

Return the KSP if the PC is Type.KSP.

getMGCoarseSolve()

Return the KSP used on the coarse grid.

getMGInterpolation(level)

Return the interpolation operator for the given level.

getMGLevels()

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.

getMGType()

Return the form of multigrid.

getOperators()

Return the matrices associated with a linear system.

getOptionsPrefix()

Return the prefix used for all the PC options.

getPythonContext()

Return the instance of the class implementing the required Python methods.

getPythonType()

Return the fully qualified Python name of the class used by the preconditioner.

getType()

Return the preconditioner type.

getUseAmat()

Return the flag to indicate if PC is applied to A or P.

matApply(x, y)

Apply the PC to many vectors stored as Mat.Type.DENSE.

reset()

Reset the PC, removing any allocated vectors and matrices.

setASMLocalSubdomains(nsd[, is_, 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_, 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.

setBDDCDirichletBoundaries(bndr)

Set the IS defining Dirichlet boundaries for the global problem.

setBDDCDirichletBoundariesLocal(bndr)

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.

setBDDCNeumannBoundaries(bndr)

Set the IS defining Neumann boundaries for the global problem.

setBDDCNeumannBoundariesLocal(bndr)

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.

setDeflationCoarseMat(mat)

Set the coarse problem matrix.

setDeflationCorrectionFactor(fact)

Set the coarse problem correction factor.

setDeflationInitOnly(flg)

Set to only perform the initialization.

setDeflationLevels(levels)

Set the maximum level of deflation nesting.

setDeflationProjectionNullSpaceMat(mat)

Set the projection null space matrix.

setDeflationReductionFactor(red)

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.

setFactorSetUpSolverType()

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.

setFromOptions()

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.

setHPDDMHasNeumannMat(has)

Set to indicate that the Mat passed to the PC is the local Neumann matrix.

setHPDDMRHSMat(B)

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.

setHYPREDiscreteGradient(mat)

Set the discrete gradient matrix.

setHYPRESetAlphaPoissonMatrix(mat)

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.

setPatchCellNumbering(sec)

Source code at petsc4py/PETSc/PC.pyx:2306

setPatchComputeFunction(function[, args, kargs])

Source code at petsc4py/PETSc/PC.pyx:2363

setPatchComputeFunctionInteriorFacets(function)

Source code at petsc4py/PETSc/PC.pyx:2370

setPatchComputeOperator(operator[, args, kargs])

Source code at petsc4py/PETSc/PC.pyx:2349

setPatchComputeOperatorInteriorFacets(operator)

Source code at petsc4py/PETSc/PC.pyx:2356

setPatchConstructType(typ[, operator, args, ...])

Source code at petsc4py/PETSc/PC.pyx:2377

setPatchDiscretisationInfo(dms, bs, ...)

Source code at petsc4py/PETSc/PC.pyx:2309

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.

setSPAIBlockSize(n)

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.

setUpOnBlocks()

Set up the PC for each block.

setUseAmat(flag)

Set to indicate to apply PC to A and not P.

view([viewer])

View the PC object.

Methods Documentation

addCompositePCType(pc_type)#

Add a PC of the given type to the composite PC.

Collective.

Parameters:

pc_type (Type | str) – The type of the preconditioner to add.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1631

appendOptionsPrefix(prefix)#

Append to the prefix used for all the PC options.

Logically collective.

Parameters:

prefix (str) – The prefix to append to the current prefix.

Return type:

None

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

apply(x, y)#

Apply the PC to a vector.

Collective.

Parameters:
  • x (Vec) – The input vector.

  • y (Vec) – The output vector, cannot be the same as x.

Return type:

None

See also

PCApply

Source code at petsc4py/PETSc/PC.pyx:577

applySymmetricLeft(x, y)#

Apply the left part of a symmetric PC to a vector.

Collective.

Parameters:
  • x (Vec) – The input vector.

  • y (Vec) – The output vector, cannot be the same as x.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:637

applySymmetricRight(x, y)#

Apply the right part of a symmetric PC to a vector.

Collective.

Parameters:
  • x (Vec) – The input vector.

  • y (Vec) – The output vector, cannot be the same as x.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:656

applyTranspose(x, y)#

Apply the transpose of the PC to a vector.

Collective.

For complex numbers this applies the non-Hermitian transpose.

Parameters:
  • x (Vec) – The input vector.

  • y (Vec) – The output vector, cannot be the same as x.

Return type:

None

See also

PCApply

Source code at petsc4py/PETSc/PC.pyx:615

create(comm=None)#

Create an empty PC.

Collective.

The default preconditioner for sparse matrices is ILU or ICC with 0 fill on one process and block Jacobi (BJACOBI) with ILU or ICC in parallel. For dense matrices it is always None.

Parameters:

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

Return type:

Self

See also

destroy, PCCreate

Source code at petsc4py/PETSc/PC.pyx:248

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:

Self

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

destroy()#

Destroy the PC that was created with create.

Collective.

See also

PCDestroy

Source code at petsc4py/PETSc/PC.pyx:234

Return type:

Self

getASMSubKSP()#

Return the local KSP object for all blocks on this process.

Not collective.

See also

PCASMGetSubKSP

Source code at petsc4py/PETSc/PC.pyx:942

Return type:

list[KSP]

getCompositePC(n)#

Return a component of the composite PC.

Not collective.

Parameters:

n (int) – The index of the PC in the composition.

Return type:

None

See also

PCCompositeGetPC

Source code at petsc4py/PETSc/PC.pyx:1610

getDM()#

Return the DM associated with the PC.

Not collective.

See also

PCGetDM

Source code at petsc4py/PETSc/PC.pyx:677

Return type:

DM

getDeflationCoarseKSP()#

Return the coarse problem KSP.

Not collective.

Source code at petsc4py/PETSc/PC.pyx:2784

Return type:

KSP

getDeflationPC()#

Return the additional preconditioner.

See also

PCDeflationGetPC

Source code at petsc4py/PETSc/PC.pyx:2799

Return type:

PC

getFactorMatrix()#

Return the factored matrix.

Not collective.

Source code at petsc4py/PETSc/PC.pyx:1422

Return type:

Mat

getFactorSolverType()#

Return the solver package used to perform the factorization.

Not collective.

Source code at petsc4py/PETSc/PC.pyx:1282

Return type:

str

getFailedReason()#

Return the reason the PC terminated.

Logically collective.

This is the maximum reason over all ranks in the PC communicator.

Source code at petsc4py/PETSc/PC.pyx:504

Return type:

FailedReason

getFailedReasonRank()#

Return the reason the PC terminated on this rank.

Not collective.

Different ranks may have different reasons.

Source code at petsc4py/PETSc/PC.pyx:521

Return type:

FailedReason

getFieldSplitSchurGetSubKSP()#

Return the KSP for the Schur complement based splits.

Source code at petsc4py/PETSc/PC.pyx:1529

Return type:

list[KSP]

getFieldSplitSubKSP()#

Return the KSP for all splits.

Source code at petsc4py/PETSc/PC.pyx:1511

Return type:

list[KSP]

getHPDDMCoarseCorrectionType()#

Return the coarse correction type.

Source code at petsc4py/PETSc/PC.pyx:2458

Return type:

HPDDMCoarseCorrectionType

getHPDDMSTShareSubKSP()#

Return true if the KSP in SLEPc ST and the subdomain solver is shared.

Source code at petsc4py/PETSc/PC.pyx:2470

Return type:

bool

getHYPREType()#

Return the Type.HYPRE type.

See also

PCHYPREGetType

Source code at petsc4py/PETSc/PC.pyx:1072

Return type:

str

getKSP()#

Return the KSP if the PC is Type.KSP.

Not collective.

See also

PCKSPGetKSP

Source code at petsc4py/PETSc/PC.pyx:1652

getMGCoarseSolve()#

Return the KSP used on the coarse grid.

Not collective.

Source code at petsc4py/PETSc/PC.pyx:1726

Return type:

KSP

getMGInterpolation(level)#

Return the interpolation operator for the given level.

Logically collective.

Parameters:

level (int) – The level where interpolation is defined from level-1 to level.

Return type:

Mat

Source code at petsc4py/PETSc/PC.pyx:1761

getMGLevels()#

Return the number of MG levels.

Not collective.

See also

PCMGGetLevels

Source code at petsc4py/PETSc/PC.pyx:1696

Return type:

int

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 to level-1.

Return type:

Vec

See also

PCMGGetRScale

Source code at petsc4py/PETSc/PC.pyx:1843

getMGRestriction(level)#

Return the restriction operator for the given level.

Logically collective.

Parameters:

level (int) – The level where restriction is defined from level to level-1.

Return type:

Mat

Source code at petsc4py/PETSc/PC.pyx:1802

getMGSmoother(level)#

Return the KSP to be used as a smoother.

Not collective.

Parameters:

level (int) – The level of the smoother.

Return type:

KSP

Source code at petsc4py/PETSc/PC.pyx:1864

getMGSmootherDown(level)#

Return the KSP to be used as a smoother before coarse grid correction.

Not collective.

Parameters:

level (int) – The level of the smoother.

Return type:

KSP

Source code at petsc4py/PETSc/PC.pyx:1885

getMGSmootherUp(level)#

Return the KSP to be used as a smoother after coarse grid correction.

Not collective.

Parameters:

level (int) – The level of the smoother.

Return type:

KSP

Source code at petsc4py/PETSc/PC.pyx:1906

getMGType()#

Return the form of multigrid.

Logically collective.

See also

PCMGGetType

Source code at petsc4py/PETSc/PC.pyx:1669

Return type:

MGType

getOperators()#

Return the matrices associated with a linear system.

Not collective.

Source code at petsc4py/PETSc/PC.pyx:400

Return type:

tuple[Mat, Mat]

getOptionsPrefix()#

Return the prefix used for all the PC options.

Not collective.

Source code at petsc4py/PETSc/PC.pyx:325

Return type:

str

getPythonContext()#

Return the instance of the class implementing the required Python methods.

Not collective.

Source code at petsc4py/PETSc/PC.pyx:775

Return type:

Any

getPythonType()#

Return the fully qualified Python name of the class used by the preconditioner.

Not collective.

Source code at petsc4py/PETSc/PC.pyx:804

Return type:

str

getType()#

Return the preconditioner type.

Not collective.

See also

setType, PCGetType

Source code at petsc4py/PETSc/PC.pyx:292

Return type:

str

getUseAmat()#

Return the flag to indicate if PC is applied to A or P.

Logically collective.

Returns:

flag – True if A is used and False if P.

Return type:

bool

Source code at petsc4py/PETSc/PC.pyx:442

matApply(x, y)#

Apply the PC to many vectors stored as Mat.Type.DENSE.

Collective.

Parameters:
  • x (Mat) – The input matrix.

  • y (Mat) – The output matrix, cannot be the same as x.

Return type:

None

See also

PCMatApply, PCApply

Source code at petsc4py/PETSc/PC.pyx:596

reset()#

Reset the PC, removing any allocated vectors and matrices.

Collective.

See also

PCReset

Source code at petsc4py/PETSc/PC.pyx:549

Return type:

None

setASMLocalSubdomains(nsd, is_=None, is_local=None)#

Set the local subdomains.

Collective.

Parameters:
  • nsd (int) – The number of subdomains for this process.

  • is_ (Sequence[IS] | None) – Defines the subdomains for this process or None to determine internally.

  • is_local (Sequence[IS] | None) – Defines the local part of the subdomains for this process, only used for PC.ASMType.RESTRICT.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:856

setASMOverlap(overlap)#

Set the overlap between a pair of subdomains.

Logically collective.

Parameters:

overlap (int) – The amount of overlap between subdomains.

Return type:

None

See also

PCASMSetOverlap

Source code at petsc4py/PETSc/PC.pyx:838

setASMSortIndices(dosort)#

Set to sort subdomain indices.

Logically collective.

Parameters:

dosort (bool) – Set to True to sort indices

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:957

setASMTotalSubdomains(nsd, is_=None, is_local=None)#

Set the subdomains for all processes.

Collective.

Parameters:
  • nsd (int) – The number of subdomains for all processes.

  • is_ (Sequence[IS] | None) – Defines the subdomains for all processes or None to determine internally.

  • is_local (Sequence[IS] | None) – Defines the local part of the subdomains for this process, only used for PC.ASMType.RESTRICT.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:899

setASMType(asmtype)#

Set the type of restriction and interpolation.

Logically collective.

Parameters:

asmtype (ASMType) – The type of ASM you wish to use.

Return type:

None

See also

PCASMSetType

Source code at petsc4py/PETSc/PC.pyx:820

setBDDCChangeOfBasisMat(T, interior=False)#

Set a user defined change of basis for degrees of freedom.

Collective.

Parameters:
  • T (Mat) – The matrix representing the change of basis.

  • interior (bool) – Enable to indicate the change of basis affects interior degrees of freedom.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2097

setBDDCCoarseningRatio(cratio)#

Set the coarsening ratio used in the multilevel version.

Logically collective.

Parameters:

cratio (int) – The coarsening ratio at the coarse level

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2152

setBDDCDirichletBoundaries(bndr)#

Set the IS defining Dirichlet boundaries for the global problem.

Collective.

Parameters:

bndr (IS) – The parallel IS defining Dirichlet boundaries.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2188

setBDDCDirichletBoundariesLocal(bndr)#

Set the IS defining Dirichlet boundaries in local ordering.

Collective.

Parameters:

bndr (IS) – The parallel IS defining Dirichlet boundaries in local ordering.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2205

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:

None

Source code at petsc4py/PETSc/PC.pyx:2061

setBDDCDivergenceMat(div, trans=False, l2l=None)#

Set the linear operator representing ∫ div(u)•p dx.

Collective.

Parameters:
  • div (Mat) – The matrix in Mat.Type.IS format.

  • trans (bool) – If True, the pressure/velocity is in the trial/test space respectively. If False the pressure/velocity is in the test/trial space.

  • l2l (IS | None) – Optional IS describing the local to local map for velocities.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2035

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:

None

Source code at petsc4py/PETSc/PC.pyx:2256

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.

Source code at petsc4py/PETSc/PC.pyx:2279

setBDDCLevels(levels)#

Set the maximum number of additional levels allowed.

Logically collective.

Parameters:

levels (int) – The maximum number of levels.

Return type:

None

See also

PCBDDCSetLevels

Source code at petsc4py/PETSc/PC.pyx:2170

setBDDCNeumannBoundaries(bndr)#

Set the IS defining Neumann boundaries for the global problem.

Collective.

Parameters:

bndr (IS) – The parallel IS defining Neumann boundaries.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2222

setBDDCNeumannBoundariesLocal(bndr)#

Set the IS defining Neumann boundaries in local ordering.

Collective.

Parameters:

bndr (IS) – The parallel IS defining Neumann boundaries in local ordering.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2239

setBDDCPrimalVerticesIS(primv)#

Set additional user defined primal vertices.

Collective.

Parameters:

primv (IS) – The IS of primal vertices in global numbering.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2118

setBDDCPrimalVerticesLocalIS(primv)#

Set additional user defined primal vertices.

Collective.

Parameters:

primv (IS) – The IS of primal vertices in local numbering.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2135

setCompositeType(ctype)#

Set the type of composite preconditioner.

Logically collective.

Parameters:

ctype (CompositeType) – The type of composition.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1592

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:

None

See also

PCSetCoordinates

Source code at petsc4py/PETSc/PC.pyx:711

setDM(dm)#

Set the DM that may be used by some preconditioners.

Logically collective.

Parameters:

dm (DM) – The DM object.

Return type:

None

See also

PCSetDM

Source code at petsc4py/PETSc/PC.pyx:694

setDeflationCoarseMat(mat)#

Set the coarse problem matrix.

Collective.

Parameters:

mat (Mat) – The coarse problem matrix.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2767

setDeflationCorrectionFactor(fact)#

Set the coarse problem correction factor.

Logically collective.

Parameters:

fact (float) – The correction factor.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2690

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.

Parameters:

flg (bool) – Enable to only initialize the preconditioner.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2632

setDeflationLevels(levels)#

Set the maximum level of deflation nesting.

Logically collective.

Parameters:

levels (int) – The maximum deflation level.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2654

setDeflationProjectionNullSpaceMat(mat)#

Set the projection null space matrix.

Collective.

Parameters:

mat (Mat) – The projection null space matrix.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2750

setDeflationReductionFactor(red)#

Set the reduction factor for the preconditioner.

Logically collective.

Parameters:

red (int) – The reduction factor or DEFAULT.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2672

setDeflationSpace(W, transpose)#

Set the deflation space matrix or its (Hermitian) transpose.

Logically collective.

Parameters:
  • W (Mat) – The deflation matrix.

  • transpose (bool) – Enable to indicate that W is an explicit transpose of the deflation matrix.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2729

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:

None

Source code at petsc4py/PETSc/PC.pyx:2708

setFactorLevels(levels)#

Set the number of levels of fill.

Logically collective.

Parameters:

levels (int) – The number of levels to fill.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1404

setFactorOrdering(ord_type=None, nzdiag=None, reuse=None)#

Set options for the matrix factorization reordering.

Logically collective.

Parameters:
  • ord_type (str | None) – The name of the matrix ordering or None to leave unchanged.

  • nzdiag (float | None) – Threshold to consider diagonal entries in the matrix as zero.

  • reuse (bool | None) – Enable to reuse the ordering of a factored matrix.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1310

setFactorPivot(zeropivot=None, inblocks=None)#

Set options for matrix factorization pivoting.

Logically collective.

Parameters:
  • zeropivot (float | None) – The size at which smaller pivots are treated as zero.

  • inblocks (bool | None) – Enable to allow pivoting while factoring in blocks.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1347

setFactorSetUpSolverType()#

Set up the factorization solver.

This can be called after KSP.setOperators or PC.setOperators, causes MatGetFactor to be called so then one may set the options for that particular factorization object.

Source code at petsc4py/PETSc/PC.pyx:1296

Return type:

None

setFactorShift(shift_type=None, amount=None)#

Set options for shifting diagonal entries of a matrix.

Parameters:
Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1376

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:

None

Source code at petsc4py/PETSc/PC.pyx:1263

setFailedReason(reason)#

Set the reason the PC terminated.

Logically collective.

Parameters:

reason (FailedReason | str) – the reason the PC terminated

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:486

setFieldSplitFields(bsize, *fields)#

Sets the elements for the field split.

Parameters:
  • bsize (int) – The block size

  • fields (Tuple[str, Sequence[int]]) – A sequence of tuples containing the split name and a sequence of integers that define the elements in the split.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1483

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:

None

Source code at petsc4py/PETSc/PC.pyx:1457

setFieldSplitSchurFactType(ctype)#

Set the type of approximate block factorization.

Collective.

Parameters:

ctype (FieldSplitSchurFactType) – The type indicating which blocks to retain.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1547

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:

None

Source code at petsc4py/PETSc/PC.pyx:1565

setFieldSplitType(ctype)#

Set the type of composition of a field split preconditioner.

Collective.

Parameters:

ctype (CompositeType) – The type of composition.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1439

setFromOptions()#

Set various PC parameters from user options.

Collective.

Source code at petsc4py/PETSc/PC.pyx:358

Return type:

None

setGAMGLevels(levels)#

Set the maximum number of levels.

Not collective.

Parameters:

levels (int) – The maximum number of levels to use.

Return type:

None

See also

PCGAMGSetNlevels

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

setGAMGSmooths(smooths)#

Set the number of smoothing steps used on all levels.

Logically collective.

Parameters:

smooths (int) – The maximum number of smooths.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1052

setGAMGType(gamgtype)#

Set the type of algorithm.

Collective.

Parameters:

gamgtype (GAMGType | str) – The type of GAMG

Return type:

None

See also

PCGAMGSetType

Source code at petsc4py/PETSc/PC.pyx:1015

setGASMOverlap(overlap)#

Set the overlap between a pair of subdomains.

Logically collective.

Parameters:

overlap (int) – The amount of overlap between subdomains.

Return type:

None

See also

PCGASMSetOverlap

Source code at petsc4py/PETSc/PC.pyx:995

setGASMType(gasmtype)#

Set the type of restriction and interpolation.

Logically collective.

Parameters:

gasmtype (GASMType) – The type of GASM.

Return type:

None

See also

PCGASMSetType

Source code at petsc4py/PETSc/PC.pyx:977

setHPDDMAuxiliaryMat(uis, uaux)#

Set the auxiliary matrix used by the preconditioner.

Parameters:
  • uis (IS) – The IS of the local auxiliary matrix

  • uaux (Mat) – The auxiliary sequential matrix

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2392

setHPDDMCoarseCorrectionType(correction_type)#

Set the coarse correction type.

Collective.

Parameters:

correction_type (HPDDMCoarseCorrectionType) – The type of coarse correction to apply.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2440

setHPDDMDeflationMat(uis, U)#

Set the deflation space used to assemble a coarse operator.

Parameters:
  • uis (IS) – The IS of the local deflation matrix.

  • U (Mat) – The deflation sequential matrix of type Mat.Type.DENSE.

Source code at petsc4py/PETSc/PC.pyx:2482

setHPDDMHasNeumannMat(has)#

Set to indicate that the Mat passed to the PC is the local Neumann matrix.

Parameters:

has (bool) – Enable to indicate the matrix is the local Neumann matrix.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2424

setHPDDMRHSMat(B)#

Set the right-hand side matrix of the preconditioner.

Parameters:

B (Mat) – The right-hand side sequential matrix.

Return type:

None

See also

PCHPDDMSetRHSMat

Source code at petsc4py/PETSc/PC.pyx:2409

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:

None

Source code at petsc4py/PETSc/PC.pyx:1244

setHYPREDiscreteCurl(mat)#

Set the discrete curl matrix.

Collective.

Parameters:

mat (Mat) – The discrete curl.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1102

setHYPREDiscreteGradient(mat)#

Set the discrete gradient matrix.

Collective.

Parameters:

mat (Mat) – The discrete gradient.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1119

setHYPRESetAlphaPoissonMatrix(mat)#

Set the vector Poisson matrix.

Collective.

Parameters:

mat (Mat) – The vector Poisson matrix.

Return type:

None

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

setHYPRESetBetaPoissonMatrix(mat=None)#

Set the Posson matrix.

Collective.

Parameters:

mat (Mat | None) – The Poisson matrix or None to turn off.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1153

setHYPRESetEdgeConstantVectors(ozz, zoz, zzo=None)#

Set the representation of the constant vector fields in the edge element basis.

Collective.

Parameters:
  • ozz (Vec) – A vector representing [1,0,0] or [1,0] in 2D.

  • zoz (Vec) – A vector representing [0,1,0] or [0,1] in 2D.

  • zzo (Vec | None) – A vector representing [0,0,1] or None in 2D.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1220

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:

None

Source code at petsc4py/PETSc/PC.pyx:1172

setHYPREType(hypretype)#

Set the Type.HYPRE type.

Parameters:

hypretype (str) – The name of the type, one of "euclid", "pilut", "parasails", "boomeramg", "ams", "ads"

See also

PCHYPRESetType

Source code at petsc4py/PETSc/PC.pyx:1084

setMGCycleType(cycle_type)#

Set the type of cycles.

Parameters:

cycle_type (MGCycleType) – The type of multigrid cycles to use.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1927

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:

None

Source code at petsc4py/PETSc/PC.pyx:1943

setMGInterpolation(level, mat)#

Set the interpolation operator for the given level.

Logically collective.

Parameters:
  • level – The level where interpolation is defined from level-1 to level.

  • mat (Mat) – The interpolation operator

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1741

setMGLevels(levels)#

Set the number of MG levels.

Parameters:

levels (int) – The number of levels

Return type:

None

See also

PCMGSetLevels

Source code at petsc4py/PETSc/PC.pyx:1710

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:
  • level (int) – The level on which to set the residual.

  • r (Vec) – The vector where the residual is stored.

Return type:

None

See also

PCMGSetR

Source code at petsc4py/PETSc/PC.pyx:2010

setMGRScale(level, rscale)#

Set the pointwise scaling for the restriction operator on the given level.

Logically collective.

Parameters:
  • level (int) – The level where restriction is defined from level to level-1.

  • rscale (Vec) – The scaling vector.

Return type:

None

See also

PCMGSetRScale

Source code at petsc4py/PETSc/PC.pyx:1823

setMGRestriction(level, mat)#

Set the restriction operator for the given level.

Logically collective.

Parameters:
  • level (int) – The level where restriction is defined from level to level-1.

  • mat (Mat) – The restriction operator

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:1782

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:
  • level (int) – The level on which to set the right-hand side.

  • rhs (Vec) – The vector where the right-hand side is stored.

Return type:

None

See also

PCMGSetRhs

Source code at petsc4py/PETSc/PC.pyx:1964

setMGType(mgtype)#

Set the form of multigrid.

Logically collective.

See also

PCMGSetType

Source code at petsc4py/PETSc/PC.pyx:1683

Parameters:

mgtype (MGType) –

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:
  • level (int) – The level on which to set the solution.

  • x (Vec) – The vector where the solution is stored.

Return type:

None

See also

PCMGSetX

Source code at petsc4py/PETSc/PC.pyx:1987

setOperators(A=None, P=None)#

Set the matrices associated with the linear system.

Logically collective.

Passing None for A or P removes the matrix that is currently used. PETSc does not reset the matrix entries of either A or P to zero after a linear solve; the user is completely responsible for matrix assembly. See Mat.zeroEntries to zero all elements of a matrix.

Parameters:
  • A (Mat | None) – the matrix which defines the linear system

  • P (Mat | None) – the matrix to be used in constructing the preconditioner, usually the same as A

Return type:

None

See also

PCSetOperators

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

setOptionsPrefix(prefix)#

Set the prefix used for all the PC options.

Logically collective.

Parameters:

prefix (str) – The prefix to prepend to all option names.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:306

setPatchCellNumbering(sec)#

Source code at petsc4py/PETSc/PC.pyx:2306

Parameters:

sec (Section) –

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

Source code at petsc4py/PETSc/PC.pyx:2363

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

Source code at petsc4py/PETSc/PC.pyx:2370

setPatchComputeOperator(operator, args=None, kargs=None)#

Source code at petsc4py/PETSc/PC.pyx:2349

setPatchComputeOperatorInteriorFacets(operator, args=None, kargs=None)#

Source code at petsc4py/PETSc/PC.pyx:2356

setPatchConstructType(typ, operator=None, args=None, kargs=None)#

Source code at petsc4py/PETSc/PC.pyx:2377

setPatchDiscretisationInfo(dms, bs, cellNodeMaps, subspaceOffsets, ghostBcNodes, globalBcNodes)#

Source code at petsc4py/PETSc/PC.pyx:2309

setPythonContext(context)#

Set the instance of the class implementing the required Python methods.

Not collective.

Source code at petsc4py/PETSc/PC.pyx:763

Parameters:

context (Any) –

Return type:

None

setPythonType(py_type)#

Set the fully qualified Python name of the class to be used.

Collective.

Source code at petsc4py/PETSc/PC.pyx:790

Parameters:

py_type (str) –

Return type:

None

setReusePreconditioner(flag)#

Set to indicate the preconditioner is to be reused.

Logically collective.

Normally if the A matrix inside a PC changes, the PC 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 and False to recompute on changes to the matrix.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:461

setSPAIBlockSize(n)#

Set the block size of the preconditioner.

Parameters:

n (int) – The block size, defaults to 1.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2566

setSPAICacheSize(size)#

Set the cache size.

Parameters:

size (int) – The size of the cache, defaults to 5.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:2582

setSPAIEpsilon(val)#

Set the tolerance for the preconditioner.

Parameters:

val (float) – The tolerance, defaults to 0.4.

Return type:

None

See also

PCSPAISetEpsilon

Source code at petsc4py/PETSc/PC.pyx:2501

setSPAIMax(maxval)#

Set the size of working buffers in the preconditioner.

Parameters:

maxval (int) – Number of entries in the work arrays to be allocated, defaults to 5000.

Return type:

None

See also

PCSPAISetMax

Source code at petsc4py/PETSc/PC.pyx:2533

setSPAIMaxNew(maxval)#

Set the maximum number of new non-zero candidates per step.

Parameters:

maxval (int) – Number of entries allowed, defaults to 5.

Return type:

None

See also

PCSPAISetMaxNew

Source code at petsc4py/PETSc/PC.pyx:2550

setSPAINBSteps(nbsteps)#

Set the maximum number of improvement steps per row.

Parameters:

nbsteps (int) – The number of steps, defaults to 5.

Return type:

None

See also

PCSPAISetNBSteps

Source code at petsc4py/PETSc/PC.pyx:2517

setSPAISp(sym)#

Set to specify a symmetric sparsity pattern.

Parameters:

sym (int) – Enable to indicate the matrix is symmetric.

Return type:

None

See also

PCSPAISetSp

Source code at petsc4py/PETSc/PC.pyx:2614

setSPAIVerbose(level)#

Set the verbosity level.

Parameters:

level (int) – The level of verbosity, defaults to 1.

Return type:

None

See also

PCSPAISetVerbose

Source code at petsc4py/PETSc/PC.pyx:2598

setType(pc_type)#

Set the preconditioner type.

Collective.

Parameters:

pc_type (Type | str) – The preconditioner type.

Return type:

None

Source code at petsc4py/PETSc/PC.pyx:273

setUp()#

Set up the internal data structures for the PC.

Collective.

See also

PCSetUp

Source code at petsc4py/PETSc/PC.pyx:537

Return type:

None

setUpOnBlocks()#

Set up the PC for each block.

Collective.

For nested preconditioners such as BJACOBI, setUp is not called on each sub-KSP when setUp is called on the outer PC. This routine ensures it is called.

Source code at petsc4py/PETSc/PC.pyx:561

Return type:

None

setUseAmat(flag)#

Set to indicate to apply PC to A and not P.

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 to TS.setRHSJacobian, TS.setIJacobian, SNES.setJacobian, KSP.setOperators or PC.setOperators not the P.

Parameters:

flag (bool) – Set True to use A and False to use P.

Return type:

None

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

view(viewer=None)#

View the PC object.

Collective.

Parameters:

viewer (Viewer | None) – The visualization context.

Return type:

None

See also

PCView

Source code at petsc4py/PETSc/PC.pyx:215