petsc4py.PETSc.KSP#

class petsc4py.PETSc.KSP#

Bases: Object

Abstract PETSc object that manages all Krylov methods.

This is the object that manages the linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators).

Notes

When a direct solver is used, but no Krylov solver is used, the KSP object is still used but with a Type.PREONLY, meaning that only application of the preconditioner is used as the linear solver.

Enumerations

ConvergedReason

KSP Converged Reason.

HPDDMType

The HPDDM Krylov solver type.

NormType

KSP norm type.

Type

KSP Type.

Methods Summary

addConvergenceTest(converged[, args, kargs, ...])

Add the function to be used to determine convergence.

appendOptionsPrefix(prefix)

Append to prefix used for all KSP options in the database.

buildResidual([r])

Return the residual of the linear system.

buildSolution([x])

Return the solution vector.

callConvergenceTest(its, rnorm)

Call the convergence test callback.

computeEigenvalues()

Compute the extreme eigenvalues for the preconditioned operator.

computeExtremeSingularValues()

Compute the extreme singular values for the preconditioned operator.

create([comm])

Create the KSP context.

createPython([context, comm])

Create a linear solver of Python type.

destroy()

Destroy KSP context.

getAppCtx()

Return the user-defined context for the linear solver.

getComputeEigenvalues()

Return flag indicating whether eigenvalues will be calculated.

getComputeSingularValues()

Return flag indicating whether singular values will be calculated.

getConvergedReason()

Use reason property.

getConvergenceHistory()

Return array containing the residual history.

getConvergenceTest()

Return the function to be used to determine convergence.

getDM()

Return the DM that may be used by some preconditioners.

getErrorIfNotConverged()

Return the flag indicating the solver will error if divergent.

getHPDDMType()

Return the Krylov solver type.

getInitialGuessKnoll()

Determine whether the KSP solver is using the Knoll trick.

getInitialGuessNonzero()

Determine whether the KSP solver uses a zero initial guess.

getIterationNumber()

Use its property.

getMonitor()

Return function used to monitor the residual.

getNormType()

Return the norm that is used for convergence testing.

getOperators()

Return the matrix associated with the linear system.

getOptionsPrefix()

Return the prefix used for all KSP options in the database.

getPC()

Return the preconditioner.

getPCSide()

Return the preconditioning side.

getPythonContext()

Return the instance of the class implementing Python methods.

getPythonType()

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

getResidualNorm()

Use norm property.

getRhs()

Return the right-hand side vector for the linear system.

getSolution()

Return the solution for the linear system to be solved.

getTolerances()

Return various tolerances used by the KSP convergence tests.

getType()

Return the KSP type as a string from the KSP object.

getWorkVecs([right, left])

Create working vectors.

logConvergenceHistory(rnorm)

Add residual to convergence history.

matSolve(B, X)

Solve a linear system with multiple right-hand sides.

matSolveTranspose(B, X)

Solve the transpose of a linear system with multiple RHS.

monitor(its, rnorm)

Run the user provided monitor routines, if they exist.

monitorCancel()

Clear all monitors for a KSP object.

reset()

Resets a KSP context.

setAppCtx(appctx)

Set the optional user-defined context for the linear solver.

setComputeEigenvalues(flag)

Set a flag to compute eigenvalues.

setComputeOperators(operators[, args, kargs])

Set routine to compute the linear operators.

setComputeRHS(rhs[, args, kargs])

Set routine to compute the right-hand side of the linear system.

setComputeSingularValues(flag)

Set flag to calculate singular values.

setConvergedReason(reason)

Use reason property.

setConvergenceHistory([length, reset])

Set the array used to hold the residual history.

setConvergenceTest(converged[, args, kargs])

Set the function to be used to determine convergence.

setDM(dm)

Set the DM that may be used by some preconditioners.

setDMActive(flag)

DM should be used to generate system matrix & RHS vector.

setErrorIfNotConverged(flag)

Cause solve to generate an error if not converged.

setFromOptions()

Set KSP options from the options database.

setGMRESRestart(restart)

Set number of iterations at which KSP restarts.

setHPDDMType(hpddm_type)

Set the Krylov solver type.

setInitialGuessKnoll(flag)

Tell solver to use PC.apply to compute the initial guess.

setInitialGuessNonzero(flag)

Tell the iterative solver that the initial guess is nonzero.

setIterationNumber(its)

Use its property.

setMonitor(monitor[, args, kargs])

Set additional function to monitor the residual.

setNormType(normtype)

Set the norm that is used for convergence testing.

setOperators([A, P])

Set matrix associated with the linear system.

setOptionsPrefix(prefix)

Set the prefix used for all KSP options in the database.

setPC(pc)

Set the preconditioner.

setPCSide(side)

Set the preconditioning side.

setPostSolve(postsolve[, args, kargs])

Set the function that is called at the end of each KSP.solve.

setPreSolve(presolve[, args, kargs])

Set the function that is called at the beginning of each KSP.solve.

setPythonContext([context])

Set the instance of the class implementing Python methods.

setPythonType(py_type)

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

setResidualNorm(rnorm)

Use norm property.

setTolerances([rtol, atol, divtol, max_it])

Set various tolerances used by the KSP convergence testers.

setType(ksp_type)

Build the KSP data structure for a particular Type.

setUp()

Set up internal data structures for an iterative solver.

setUpOnBlocks()

Set up the preconditioner for each block in a block method.

setUseFischerGuess(model, size)

Use the Paul Fischer algorithm to compute initial guesses.

solve(b, x)

Solve the linear system.

solveTranspose(b, x)

Solve the transpose of a linear system.

view([viewer])

Print the KSP data structure.

Attributes Summary

appctx

The solver application context.

atol

The absolute tolerance of the solver.

divtol

The divergence tolerance of the solver.

dm

The solver DM.

guess_knoll

Whether solver uses Knoll trick.

guess_nonzero

Whether guess is non-zero.

history

The convergence history of the solver.

is_converged

Boolean indicating if the solver has converged.

is_diverged

Boolean indicating if the solver has failed.

is_iterating

Boolean indicating if the solver has not converged yet.

its

The current number of iterations the solver has taken.

mat_op

The system matrix operator.

mat_pc

The preconditioner operator.

max_it

The maximum number of iteration the solver may take.

norm

The norm of the residual at the current iteration.

norm_type

The norm used by the solver.

pc

The PC of the solver.

pc_side

The side on which preconditioning is performed.

reason

The converged reason.

rtol

The relative tolerance of the solver.

vec_rhs

The right-hand side vector.

vec_sol

The solution vector.

Methods Documentation

addConvergenceTest(converged, args=None, kargs=None, prepend=False)#

Add the function to be used to determine convergence.

Logically collective.

Parameters:
  • converged (KSPConvergenceTestFunction) – Callback which computes the convergence.

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

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

  • prepend (bool) – Whether to prepend this call before the default convergence test or call it after.

Return type:

None

Notes

Cannot be mixed with a call to setConvergenceTest. It can only be called once. If called multiple times, it will generate an error.

Source code at petsc4py/PETSc/KSP.pyx:1057

appendOptionsPrefix(prefix)#

Append to prefix used for all KSP options in the database.

Logically collective.

Parameters:

prefix (str | None) – The options prefix to append.

Return type:

None

Notes

A hyphen (-) must NOT be given at the beginning of the prefix name. The first character of all runtime options is AUTOMATICALLY the hyphen.

Source code at petsc4py/PETSc/KSP.pyx:575

buildResidual(r=None)#

Return the residual of the linear system.

Collective.

Parameters:

r (Vec | None) – Optional vector to use for the result.

Return type:

Vec

Source code at petsc4py/PETSc/KSP.pyx:2041

buildSolution(x=None)#

Return the solution vector.

Collective.

Parameters:

x (Vec | None) – Optional vector to store the solution.

Return type:

Vec

Source code at petsc4py/PETSc/KSP.pyx:2019

callConvergenceTest(its, rnorm)#

Call the convergence test callback.

Collective.

Parameters:
  • its (int) – Number of iterations.

  • rnorm (float) – The residual norm.

Return type:

None

Notes

This functionality is implemented in petsc4py.

Source code at petsc4py/PETSc/KSP.pyx:1113

computeEigenvalues()#

Compute the extreme eigenvalues for the preconditioned operator.

Not collective.

Source code at petsc4py/PETSc/KSP.pyx:2063

Return type:

ArrayComplex

computeExtremeSingularValues()#

Compute the extreme singular values for the preconditioned operator.

Collective.

Returns:

  • smax (float) – The maximum singular value.

  • smin (float) – The minimum singular value.

Return type:

tuple[float, float]

Source code at petsc4py/PETSc/KSP.pyx:2086

create(comm=None)#

Create the KSP context.

Collective.

See also

KSPCreate

Source code at petsc4py/PETSc/KSP.pyx:457

Parameters:

comm (Comm | None)

Return type:

Self

createPython(context=None, comm=None)#

Create a linear solver 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/KSP.pyx:2132

destroy()#

Destroy KSP context.

Collective.

See also

KSPDestroy

Source code at petsc4py/PETSc/KSP.pyx:444

Return type:

Self

getAppCtx()#

Return the user-defined context for the linear solver.

Not collective.

See also

setAppCtx

Source code at petsc4py/PETSc/KSP.pyx:640

Return type:

Any

getComputeEigenvalues()#

Return flag indicating whether eigenvalues will be calculated.

Not collective.

Return the flag indicating that the extreme eigenvalues values will be calculated via a Lanczos or Arnoldi process as the linear system is solved.

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

Return type:

bool

getComputeSingularValues()#

Return flag indicating whether singular values will be calculated.

Not collective.

Return the flag indicating whether the extreme singular values will be calculated via a Lanczos or Arnoldi process as the linear system is solved.

Source code at petsc4py/PETSc/KSP.pyx:1466

Return type:

bool

getConvergedReason()#

Use reason property.

Source code at petsc4py/PETSc/KSP.pyx:1868

Return type:

ConvergedReason

getConvergenceHistory()#

Return array containing the residual history.

Not collective.

Source code at petsc4py/PETSc/KSP.pyx:1180

Return type:

ArrayReal

getConvergenceTest()#

Return the function to be used to determine convergence.

Logically collective.

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

Return type:

KSPConvergenceTestFunction

getDM()#

Return the DM that may be used by some preconditioners.

Not collective.

See also

PETSc.KSP, DM, KSPGetDM

Source code at petsc4py/PETSc/KSP.pyx:654

Return type:

DM

getErrorIfNotConverged()#

Return the flag indicating the solver will error if divergent.

Not collective.

Source code at petsc4py/PETSc/KSP.pyx:1924

Return type:

bool

getHPDDMType()#

Return the Krylov solver type.

Not collective.

See also

KSPHPDDMGetType

Source code at petsc4py/PETSc/KSP.pyx:1892

Return type:

HPDDMType

getInitialGuessKnoll()#

Determine whether the KSP solver is using the Knoll trick.

Not collective.

This uses the Knoll trick; using PC.apply to compute the initial guess.

Source code at petsc4py/PETSc/KSP.pyx:1542

Return type:

bool

getInitialGuessNonzero()#

Determine whether the KSP solver uses a zero initial guess.

Not collective.

Source code at petsc4py/PETSc/KSP.pyx:1508

Return type:

bool

getIterationNumber()#

Use its property.

Source code at petsc4py/PETSc/KSP.pyx:1846

Return type:

int

getMonitor()#

Return function used to monitor the residual.

Not collective.

Source code at petsc4py/PETSc/KSP.pyx:1258

Return type:

KSPMonitorFunction

getNormType()#

Return the norm that is used for convergence testing.

Not collective.

Source code at petsc4py/PETSc/KSP.pyx:1382

Return type:

NormType

getOperators()#

Return the matrix associated with the linear system.

Collective.

Return the matrix associated with the linear system and a (possibly) different one used to construct the preconditioner.

Returns:

  • A (Mat) – Matrix that defines the linear system.

  • P (Mat) – Matrix to be used in constructing the preconditioner, usually the same as A.

Return type:

tuple[Mat, Mat]

Source code at petsc4py/PETSc/KSP.pyx:850

getOptionsPrefix()#

Return the prefix used for all KSP options in the database.

Not collective.

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

Return type:

str

getPC()#

Return the preconditioner.

Not collective.

See also

PETSc.KSP, setPC, KSPGetPC

Source code at petsc4py/PETSc/KSP.pyx:897

Return type:

PC

getPCSide()#

Return the preconditioning side.

Not collective.

Source code at petsc4py/PETSc/KSP.pyx:1341

Return type:

Side

getPythonContext()#

Return the instance of the class implementing Python methods.

Not collective.

Source code at petsc4py/PETSc/KSP.pyx:2173

Return type:

Any

getPythonType()#

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

Not collective.

Source code at petsc4py/PETSc/KSP.pyx:2203

Return type:

str

getResidualNorm()#

Use norm property.

Source code at petsc4py/PETSc/KSP.pyx:1857

Return type:

float

getRhs()#

Return the right-hand side vector for the linear system.

Not collective.

See also

KSPGetRhs

Source code at petsc4py/PETSc/KSP.pyx:1938

Return type:

Vec

getSolution()#

Return the solution for the linear system to be solved.

Not collective.

Note that this may not be the solution that is stored during the iterative process.

See also

KSPGetSolution

Source code at petsc4py/PETSc/KSP.pyx:1953

Return type:

Vec

getTolerances()#

Return various tolerances used by the KSP convergence tests.

Not collective.

Return the relative, absolute, divergence, and maximum iteration tolerances used by the default KSP convergence tests.

Returns:

  • rtol (float) – The relative convergence tolerance

  • atol (float) – The absolute convergence tolerance

  • dtol (float) – The divergence tolerance

  • maxits (int) – Maximum number of iterations

Return type:

tuple[float, float, float, int]

Source code at petsc4py/PETSc/KSP.pyx:963

getType()#

Return the KSP type as a string from the KSP object.

Not collective.

See also

KSPGetType

Source code at petsc4py/PETSc/KSP.pyx:509

Return type:

str

getWorkVecs(right=None, left=None)#

Create working vectors.

Collective.

Parameters:
  • right (int | None) – Number of right hand vectors to allocate.

  • left (int | None) – Number of left hand vectors to allocate.

Returns:

  • R (list of Vec) – List of correctly allocated right hand vectors.

  • L (list of Vec) – List of correctly allocated left hand vectors.

Return type:

tuple[list[Vec], list[Vec]] | list[Vec] | None

Source code at petsc4py/PETSc/KSP.pyx:1971

logConvergenceHistory(rnorm)#

Add residual to convergence history.

Logically collective.

Parameters:

rnorm (float) – Residual norm to be added to convergence history.

Return type:

None

Source code at petsc4py/PETSc/KSP.pyx:1195

matSolve(B, X)#

Solve a linear system with multiple right-hand sides.

Collective.

These are stored as a Mat.Type.DENSE. Unlike solve, B and X must be different matrices.

Parameters:
  • B (Mat) – Block of right-hand sides.

  • X (Mat) – Block of solutions.

Return type:

None

See also

solve, KSPMatSolve

Source code at petsc4py/PETSc/KSP.pyx:1800

matSolveTranspose(B, X)#

Solve the transpose of a linear system with multiple RHS.

Collective.

Parameters:
  • B (Mat) – Block of right-hand sides.

  • X (Mat) – Block of solutions.

Return type:

None

Source code at petsc4py/PETSc/KSP.pyx:1822

monitor(its, rnorm)#

Run the user provided monitor routines, if they exist.

Collective.

Notes

This routine is called by the KSP implementations. It does not typically need to be called by the user.

Source code at petsc4py/PETSc/KSP.pyx:1286

Parameters:
Return type:

None

monitorCancel()#

Clear all monitors for a KSP object.

Logically collective.

Source code at petsc4py/PETSc/KSP.pyx:1271

Return type:

None

reset()#

Resets a KSP context.

Collective.

Resets a KSP context to the kspsetupcalled = 0 state and removes any allocated Vecs and Mats.

See also

KSPReset

Source code at petsc4py/PETSc/KSP.pyx:1599

Return type:

None

setAppCtx(appctx)#

Set the optional user-defined context for the linear solver.

Not collective.

Parameters:

appctx (Any) – The user defined context

Return type:

None

Notes

The user context is a way for users to attach any information to the KSP that they may need later when interacting with the solver.

See also

getAppCtx

Source code at petsc4py/PETSc/KSP.pyx:617

setComputeEigenvalues(flag)#

Set a flag to compute eigenvalues.

Logically collective.

Set a flag so that the extreme eigenvalues values will be calculated via a Lanczos or Arnoldi process as the linear system is solved.

Parameters:

flag (bool) – Boolean whether to compute eigenvalues (or not).

Return type:

None

Notes

Currently this option is not valid for all iterative methods.

Source code at petsc4py/PETSc/KSP.pyx:1396

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

Set routine to compute the linear operators.

Logically collective.

Parameters:
  • operators (KSPOperatorsFunction) – Function which computes the operators.

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

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

Return type:

None

Notes

The user provided function Operators will be called automatically at the very next call to solve. It will NOT be called at future solve calls unless either setComputeOperators or setOperators is called before that solve is called. This allows the same system to be solved several times with different right-hand side functions, but is a confusing API since one might expect it to be called for each solve.

To reuse the same preconditioner for the next solve and not compute a new one based on the most recently computed matrix call KSPSetReusePreconditioner.

Source code at petsc4py/PETSc/KSP.pyx:764

setComputeRHS(rhs, args=None, kargs=None)#

Set routine to compute the right-hand side of the linear system.

Logically collective.

Parameters:
  • rhs (KSPRHSFunction) – Function which computes the right-hand side.

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

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

Return type:

None

Notes

The routine you provide will be called each time you call solve to prepare the new right-hand side for that solve.

Source code at petsc4py/PETSc/KSP.pyx:730

setComputeSingularValues(flag)#

Set flag to calculate singular values.

Logically collective.

Set a flag so that the extreme singular values will be calculated via a Lanczos or Arnoldi process as the linear system is solved.

Parameters:

flag (bool) – Boolean whether to compute singular values (or not).

Return type:

None

Notes

Currently this option is not valid for all iterative methods.

Source code at petsc4py/PETSc/KSP.pyx:1440

setConvergedReason(reason)#

Use reason property.

Source code at petsc4py/PETSc/KSP.pyx:1863

Parameters:

reason (ConvergedReason)

Return type:

None

setConvergenceHistory(length=None, reset=False)#

Set the array used to hold the residual history.

Not collective.

If set, this array will contain the residual norms computed at each iteration of the solver.

Parameters:
  • length (int | None) – Length of array to store history in.

  • reset (bool) – True indicates the history counter is reset to zero for each new linear solve.

Return type:

None

Notes

If length is not provided or None then a default array of length 10000 is allocated.

If the array is not long enough then once the iterations is longer than the array length solve stops recording the history.

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

setConvergenceTest(converged, args=None, kargs=None)#

Set the function to be used to determine convergence.

Logically collective.

Parameters:
Return type:

None

Notes

Must be called after the KSP type has been set so put this after a call to setType, or setFromOptions.

The default is a combination of relative and absolute tolerances. The residual value that is tested may be an approximation; routines that need exact values should compute them.

Source code at petsc4py/PETSc/KSP.pyx:992

setDM(dm)#

Set the DM that may be used by some preconditioners.

Logically collective.

Parameters:

dm (DM) – The DM object, cannot be None.

Return type:

None

Notes

If this is used then the KSP will attempt to use the DM to create the matrix and use the routine set with DM.setKSPComputeOperators. Use setDMActive(False) to instead use the matrix you have provided with setOperators.

A DM can only be used for solving one problem at a time because information about the problem is stored on the DM, even when not using interfaces like DM.setKSPComputeOperators. Use DM.clone to get a distinct DM when solving different problems using the same function space.

Source code at petsc4py/PETSc/KSP.pyx:671

setDMActive(flag)#

DM should be used to generate system matrix & RHS vector.

Logically collective.

Parameters:

flag (bool) – Whether to use the DM.

Return type:

None

Notes

By default setDM sets the DM as active, call setDMActive(False) after setDM(dm) to not have the KSP object use the DM to generate the matrices.

Source code at petsc4py/PETSc/KSP.pyx:704

setErrorIfNotConverged(flag)#

Cause solve to generate an error if not converged.

Logically collective.

Parameters:

flag (bool) – True enables this behavior.

Return type:

None

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

setFromOptions()#

Set KSP options from the options database.

Collective.

This routine must be called before setUp if the user is to be allowed to set the Krylov type.

Source code at petsc4py/PETSc/KSP.pyx:600

Return type:

None

setGMRESRestart(restart)#

Set number of iterations at which KSP restarts.

Logically collective.

Suitable KSPs are: KSPGMRES, KSPFGMRES and KSPLGMRES.

Parameters:

restart (int) – Integer restart value.

Return type:

None

Source code at petsc4py/PETSc/KSP.pyx:2110

setHPDDMType(hpddm_type)#

Set the Krylov solver type.

Collective.

Parameters:

hpddm_type (HPDDMType) – The type of Krylov solver to use.

Return type:

None

See also

KSPHPDDMSetType

Source code at petsc4py/PETSc/KSP.pyx:1874

setInitialGuessKnoll(flag)#

Tell solver to use PC.apply to compute the initial guess.

Logically collective.

This is the Knoll trick.

Parameters:

flag (bool) – True uses Knoll trick.

Return type:

None

Source code at petsc4py/PETSc/KSP.pyx:1522

setInitialGuessNonzero(flag)#

Tell the iterative solver that the initial guess is nonzero.

Logically collective.

Otherwise KSP assumes the initial guess is to be zero (and thus zeros it out before solving).

Parameters:

flag (bool) – True indicates the guess is non-zero, False indicates the guess is zero.

Return type:

None

Source code at petsc4py/PETSc/KSP.pyx:1486

setIterationNumber(its)#

Use its property.

Source code at petsc4py/PETSc/KSP.pyx:1841

Parameters:

its (int)

Return type:

None

setMonitor(monitor, args=None, kargs=None)#

Set additional function to monitor the residual.

Logically collective.

Set an ADDITIONAL function to be called at every iteration to monitor the residual/error etc.

Parameters:
  • monitor (KSPMonitorFunction) – Callback which monitors the convergence.

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

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

Return type:

None

Notes

The default is to do nothing. To print the residual, or preconditioned residual if setNormType(NORM_PRECONDITIONED) was called, use monitor as the monitoring routine, with a PETSc.Viewer.ASCII as the context.

Several different monitoring routines may be set by calling setMonitor multiple times; all will be called in the order in which they were set.

Source code at petsc4py/PETSc/KSP.pyx:1211

setNormType(normtype)#

Set the norm that is used for convergence testing.

Logically collective.

Parameters:

normtype (NormType) – The norm type to use (see NormType).

Return type:

None

Notes

Not all combinations of preconditioner side (see setPCSide) and norm type are supported by all Krylov methods. If only one is set, PETSc tries to automatically change the other to find a compatible pair. If no such combination is supported, PETSc will generate an error.

Source code at petsc4py/PETSc/KSP.pyx:1355

setOperators(A=None, P=None)#

Set matrix associated with the linear system.

Collective.

Set the matrix associated with the linear system and a (possibly) different one from which the preconditioner will be built.

Parameters:
  • A (Mat | None) – Matrix that defines the linear system.

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

Return type:

None

Notes

If you know the operator A has a null space you can use Mat.setNullSpace and Mat.setTransposeNullSpace to supply the null space to A and the KSP solvers will automatically use that null space as needed during the solution process.

All future calls to setOperators must use the same size matrices!

Passing None for A or P removes the matrix that is currently used.

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

setOptionsPrefix(prefix)#

Set the prefix used for all KSP options in the database.

Logically collective.

Parameters:

prefix (str | None) – The options prefix.

Return type:

None

Notes

A hyphen (-) must NOT be given at the beginning of the prefix name. The first character of all runtime options is AUTOMATICALLY the hyphen. For example, to distinguish between the runtime options for two different KSP contexts, one could call ` KSPSetOptionsPrefix(ksp1, "sys1_") KSPSetOptionsPrefix(ksp2, "sys2_") `

This would enable use of different options for each system, such as ` -sys1_ksp_type gmres -sys1_ksp_rtol 1.e-3 -sys2_ksp_type bcgs  -sys2_ksp_rtol 1.e-4 `

Source code at petsc4py/PETSc/KSP.pyx:523

setPC(pc)#

Set the preconditioner.

Collective.

Set the preconditioner to be used to calculate the application of the preconditioner on a vector.

Parameters:

pc (PC) – The preconditioner object

Return type:

None

See also

PETSc.KSP, getPC, KSPSetPC

Source code at petsc4py/PETSc/KSP.pyx:877

setPCSide(side)#

Set the preconditioning side.

Logically collective.

Parameters:

side (Side) – The preconditioning side (see PC.Side).

Return type:

None

Notes

Left preconditioning is used by default for most Krylov methods except Type.FGMRES which only supports right preconditioning.

For methods changing the side of the preconditioner changes the norm type that is used, see setNormType.

Symmetric preconditioning is currently available only for the Type.QCG method. Note, however, that symmetric preconditioning can be emulated by using either right or left preconditioning and a pre or post processing step.

Setting the PC side often affects the default norm type. See setNormType for details.

Source code at petsc4py/PETSc/KSP.pyx:1307

setPostSolve(postsolve, args=None, kargs=None)#

Set the function that is called at the end of each KSP.solve.

Logically collective.

Parameters:
  • postsolve (KSPPostSolveFunction | None) – The callback function.

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

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

Return type:

None

Source code at petsc4py/PETSc/KSP.pyx:1662

setPreSolve(presolve, args=None, kargs=None)#

Set the function that is called at the beginning of each KSP.solve.

Logically collective.

Parameters:
  • presolve (KSPPreSolveFunction | None) – The callback function.

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

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

Return type:

None

Source code at petsc4py/PETSc/KSP.pyx:1629

setPythonContext(context=None)#

Set the instance of the class implementing Python methods.

Not collective.

Source code at petsc4py/PETSc/KSP.pyx:2161

Parameters:

context (Any | None)

Return type:

None

setPythonType(py_type)#

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

Collective.

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

Parameters:

py_type (str)

Return type:

None

setResidualNorm(rnorm)#

Use norm property.

Source code at petsc4py/PETSc/KSP.pyx:1852

Parameters:

rnorm (float)

Return type:

None

setTolerances(rtol=None, atol=None, divtol=None, max_it=None)#

Set various tolerances used by the KSP convergence testers.

Logically collective.

Set the relative, absolute, divergence, and maximum iteration tolerances used by the default KSP convergence testers.

Parameters:
  • rtol (float | None) – The relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm.

  • atol (float | None) – The absolute convergence tolerance absolute size of the (possibly preconditioned) residual norm.

  • dtol – The divergence tolerance, amount (possibly preconditioned) residual norm can increase before KSPConvergedDefault concludes that the method is diverging.

  • max_it (int | None) – Maximum number of iterations to use.

  • divtol (float | None)

Return type:

None

Notes

Use None to retain the default value of any of the tolerances.

Source code at petsc4py/PETSc/KSP.pyx:914

setType(ksp_type)#

Build the KSP data structure for a particular Type.

Logically collective.

Parameters:

ksp_type (Type | str) – KSP Type object

Return type:

None

Notes

See Type for available methods (for instance, Type.CG or Type.GMRES).

Normally, it is best to use the setFromOptions command and then set the KSP type from the options database rather than by using this routine. Using the options database provides the user with maximum flexibility in evaluating the many different Krylov methods. This method is provided for those situations where it is necessary to set the iterative solver independently of the command line or options database. This might be the case, for example, when the choice of iterative solver changes during the execution of the program, and the user’s application is taking responsibility for choosing the appropriate method. In other words, this routine is not for beginners.

See also

KSPSetType

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

setUp()#

Set up internal data structures for an iterative solver.

Collective.

See also

KSPSetUp

Source code at petsc4py/PETSc/KSP.pyx:1587

Return type:

None

setUpOnBlocks()#

Set up the preconditioner for each block in a block method.

Collective.

Methods include: block Jacobi, block Gauss-Seidel, and overlapping Schwarz methods.

See also

KSPSetUpOnBlocks

Source code at petsc4py/PETSc/KSP.pyx:1614

Return type:

None

setUseFischerGuess(model, size)#

Use the Paul Fischer algorithm to compute initial guesses.

Logically collective.

Use the Paul Fischer algorithm or its variants to compute initial guesses for a set of solves with related right hand sides.

Parameters:
  • model (int) – Use model 1, model 2, model 3, any other number to turn it off.

  • size (int) – Size of subspace used to generate initial guess.

Return type:

None

Source code at petsc4py/PETSc/KSP.pyx:1559

solve(b, x)#

Solve the linear system.

Collective.

Parameters:
  • b (Vec) – Right hand side vector.

  • x (Vec) – Solution vector.

Return type:

None

Notes

If one uses setDM then x or b need not be passed. Use getSolution to access the solution in this case.

The operator is specified with setOperators.

solve will normally return without generating an error regardless of whether the linear system was solved or if constructing the preconditioner failed. Call getConvergedReason to determine if the solver converged or failed and why. The option -ksp_error_if_not_converged or function setErrorIfNotConverged will cause solve to error as soon as an error occurs in the linear solver. In inner solves, DIVERGED_MAX_IT is not treated as an error because when using nested solvers it may be fine that inner solvers in the preconditioner do not converge during the solution process.

The number of iterations can be obtained from its.

If you provide a matrix that has a Mat.setNullSpace and Mat.setTransposeNullSpace this will use that information to solve singular systems in the least squares sense with a norm minimizing solution.

Ax = b where b = bₚ + bₜ where bₜ is not in the range of A (and hence by the fundamental theorem of linear algebra is in the nullspace(Aᵀ), see Mat.setNullSpace.

KSP first removes bₜ producing the linear system Ax = bₚ (which has multiple solutions) and solves this to find the ∥x∥ minimizing solution (and hence it finds the solution x orthogonal to the nullspace(A). The algorithm is simply in each iteration of the Krylov method we remove the nullspace(A) from the search direction thus the solution which is a linear combination of the search directions has no component in the nullspace(A).

We recommend always using Type.GMRES for such singular systems. If nullspace(A) = nullspace(Aᵀ) (note symmetric matrices always satisfy this property) then both left and right preconditioning will work If nullspace(A) != nullspace(Aᵀ) then left preconditioning will work but right preconditioning may not work (or it may).

If using a direct method (e.g., via the KSP solver Type.PREONLY and a preconditioner such as PC.Type.LU or PC.Type.ILU, then its=1. See setTolerances for more details.

Understanding Convergence

The routines setMonitor and computeEigenvalues provide information on additional options to monitor convergence and print eigenvalue information.

Source code at petsc4py/PETSc/KSP.pyx:1695

solveTranspose(b, x)#

Solve the transpose of a linear system.

Collective.

Parameters:
  • b (Vec) – Right hand side vector.

  • x (Vec) – Solution vector.

Return type:

None

Notes

For complex numbers this solve the non-Hermitian transpose system.

Source code at petsc4py/PETSc/KSP.pyx:1776

view(viewer=None)#

Print the KSP data structure.

Collective.

Parameters:

viewer (Viewer | None) – Viewer used to display the KSP.

Return type:

None

See also

KSPView

Source code at petsc4py/PETSc/KSP.pyx:425

Attributes Documentation

appctx#

The solver application context.

Source code at petsc4py/PETSc/KSP.pyx:2220

atol#

The absolute tolerance of the solver.

Source code at petsc4py/PETSc/KSP.pyx:2313

divtol#

The divergence tolerance of the solver.

Source code at petsc4py/PETSc/KSP.pyx:2321

dm#

The solver DM.

Source code at petsc4py/PETSc/KSP.pyx:2230

guess_knoll#

Whether solver uses Knoll trick.

Source code at petsc4py/PETSc/KSP.pyx:2272

guess_nonzero#

Whether guess is non-zero.

Source code at petsc4py/PETSc/KSP.pyx:2264

history#

The convergence history of the solver.

Source code at petsc4py/PETSc/KSP.pyx:2355

is_converged#

Boolean indicating if the solver has converged.

Source code at petsc4py/PETSc/KSP.pyx:2375

is_diverged#

Boolean indicating if the solver has failed.

Source code at petsc4py/PETSc/KSP.pyx:2380

is_iterating#

Boolean indicating if the solver has not converged yet.

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

its#

The current number of iterations the solver has taken.

Source code at petsc4py/PETSc/KSP.pyx:2339

mat_op#

The system matrix operator.

Source code at petsc4py/PETSc/KSP.pyx:2252

mat_pc#

The preconditioner operator.

Source code at petsc4py/PETSc/KSP.pyx:2257

max_it#

The maximum number of iteration the solver may take.

Source code at petsc4py/PETSc/KSP.pyx:2329

norm#

The norm of the residual at the current iteration.

Source code at petsc4py/PETSc/KSP.pyx:2347

norm_type#

The norm used by the solver.

Source code at petsc4py/PETSc/KSP.pyx:2295

pc#

The PC of the solver.

Source code at petsc4py/PETSc/KSP.pyx:2282

pc_side#

The side on which preconditioning is performed.

Source code at petsc4py/PETSc/KSP.pyx:2287

reason#

The converged reason.

Source code at petsc4py/PETSc/KSP.pyx:2362

rtol#

The relative tolerance of the solver.

Source code at petsc4py/PETSc/KSP.pyx:2305

vec_rhs#

The right-hand side vector.

Source code at petsc4py/PETSc/KSP.pyx:2245

vec_sol#

The solution vector.

Source code at petsc4py/PETSc/KSP.pyx:2240