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
KSP Converged Reason.
The HPDDM Krylov solver type.
KSP norm 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.
Compute the extreme eigenvalues for the preconditioned operator.
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.
Return the user-defined context for the linear solver.
Return the CG objective function value.
Return flag indicating whether eigenvalues will be calculated.
Return flag indicating whether singular values will be calculated.
Use
reason
property.Return array containing the residual history.
Return the function to be used to determine convergence.
getDM
()Return the
DM
that may be used by some preconditioners.Return the flag indicating the solver will error if divergent.
Return the Krylov solver type.
Determine whether the KSP solver is using the Knoll trick.
Determine whether the KSP solver uses a zero initial guess.
Use
its
property.Return function used to monitor the residual.
Return the norm that is used for convergence testing.
Return the matrix associated with the linear system.
Return the prefix used for all
KSP
options in the database.getPC
()Return the preconditioner.
Return the preconditioning side.
Return the instance of the class implementing Python methods.
Return the fully qualified Python name of the class used by the solver.
Use
norm
property.getRhs
()Return the right-hand side vector for the linear system.
Return the solution for the linear system to be solved.
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.
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.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)setUp
()Set up internal data structures for an iterative solver.
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
The solver application context.
The absolute tolerance of the solver.
The divergence tolerance of the solver.
The solver
DM
.Whether solver uses Knoll trick.
Whether guess is non-zero.
The convergence history of the solver.
Boolean indicating if the solver has converged.
Boolean indicating if the solver has failed.
Boolean indicating if the solver has not converged yet.
The current number of iterations the solver has taken.
The system matrix operator.
The preconditioner operator.
The maximum number of iteration the solver may take.
The norm of the residual at the current iteration.
The norm used by the solver.
The
PC
of the solver.The side on which preconditioning is performed.
The converged reason.
The relative tolerance of the solver.
The right-hand side vector.
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:
Notes
Cannot be mixed with a call to
setConvergenceTest
. It can only be called once. If called multiple times, it will generate an error.
- appendOptionsPrefix(prefix)#
Append to prefix used for all
KSP
options in the database.Logically collective.
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.
See also
- buildResidual(r=None)#
Return the residual of the linear system.
Collective.
See also
- buildSolution(x=None)#
Return the solution vector.
Collective.
See also
- callConvergenceTest(its, rnorm)#
Call the convergence test callback.
Collective.
Notes
This functionality is implemented in petsc4py.
- computeEigenvalues()#
Compute the extreme eigenvalues for the preconditioned operator.
Not collective.
See also
Source code at petsc4py/PETSc/KSP.pyx:2085
- Return type:
- computeExtremeSingularValues()#
Compute the extreme singular values for the preconditioned operator.
Collective.
- Returns:
- Return type:
See also
- create(comm=None)#
Create the KSP context.
Collective.
See also
- 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:
- destroy()#
Destroy KSP context.
Collective.
See also
Source code at petsc4py/PETSc/KSP.pyx:444
- Return type:
- getAppCtx()#
Return the user-defined context for the linear solver.
Not collective.
See also
Source code at petsc4py/PETSc/KSP.pyx:640
- Return type:
- getCGObjectiveValue()#
Return the CG objective function value.
Not collective.
See also
Source code at petsc4py/PETSc/KSP.pyx:1882
- Return type:
- 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:1430
- Return type:
- 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:1474
- Return type:
- getConvergedReason()#
Use
reason
property.Source code at petsc4py/PETSc/KSP.pyx:1876
- Return type:
- getConvergenceHistory()#
Return array containing the residual history.
Not collective.
See also
Source code at petsc4py/PETSc/KSP.pyx:1188
- Return type:
- getConvergenceTest()#
Return the function to be used to determine convergence.
Logically collective.
Source code at petsc4py/PETSc/KSP.pyx:1108
- Return type:
- getDM()#
Return the
DM
that may be used by some preconditioners.Not collective.
Source code at petsc4py/PETSc/KSP.pyx:654
- Return type:
- getErrorIfNotConverged()#
Return the flag indicating the solver will error if divergent.
Not collective.
See also
Source code at petsc4py/PETSc/KSP.pyx:1946
- Return type:
- getHPDDMType()#
Return the Krylov solver type.
Not collective.
See also
Source code at petsc4py/PETSc/KSP.pyx:1914
- Return type:
- 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.See also
Source code at petsc4py/PETSc/KSP.pyx:1550
- Return type:
- getInitialGuessNonzero()#
Determine whether the KSP solver uses a zero initial guess.
Not collective.
See also
Source code at petsc4py/PETSc/KSP.pyx:1516
- Return type:
- getIterationNumber()#
Use
its
property.Source code at petsc4py/PETSc/KSP.pyx:1854
- Return type:
- getMonitor()#
Return function used to monitor the residual.
Not collective.
Source code at petsc4py/PETSc/KSP.pyx:1266
- Return type:
- getNormType()#
Return the norm that is used for convergence testing.
Not collective.
See also
Source code at petsc4py/PETSc/KSP.pyx:1390
- Return type:
- 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:
- Return type:
See also
- getOptionsPrefix()#
Return the prefix used for all
KSP
options in the database.Not collective.
See also
Source code at petsc4py/PETSc/KSP.pyx:561
- Return type:
- getPC()#
Return the preconditioner.
Not collective.
Source code at petsc4py/PETSc/KSP.pyx:897
- Return type:
- getPCSide()#
Return the preconditioning side.
Not collective.
Source code at petsc4py/PETSc/KSP.pyx:1349
- Return type:
- getPythonContext()#
Return the instance of the class implementing Python methods.
Not collective.
Source code at petsc4py/PETSc/KSP.pyx:2195
- Return type:
- getPythonType()#
Return the fully qualified Python name of the class used by the solver.
Not collective.
Source code at petsc4py/PETSc/KSP.pyx:2225
- Return type:
- getResidualNorm()#
Use
norm
property.Source code at petsc4py/PETSc/KSP.pyx:1865
- Return type:
- getRhs()#
Return the right-hand side vector for the linear system.
Not collective.
See also
Source code at petsc4py/PETSc/KSP.pyx:1960
- Return type:
- 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
Source code at petsc4py/PETSc/KSP.pyx:1975
- Return type:
- 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:
- Return type:
See also
- getType()#
Return the KSP type as a string from the
KSP
object.Not collective.
See also
Source code at petsc4py/PETSc/KSP.pyx:509
- Return type:
- getWorkVecs(right=None, left=None)#
Create working vectors.
Collective.
- Parameters:
- Returns:
- Return type:
- logConvergenceHistory(rnorm)#
Add residual to convergence history.
Logically collective.
- matSolve(B, X)#
Solve a linear system with multiple right-hand sides.
Collective.
These are stored as a
Mat.Type.DENSE
. Unlikesolve
,B
andX
must be different matrices.See also
- matSolveTranspose(B, X)#
Solve the transpose of a linear system with multiple RHS.
Collective.
See also
- 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.See also
- monitorCancel()#
Clear all monitors for a
KSP
object.Logically collective.
Source code at petsc4py/PETSc/KSP.pyx:1279
- Return type:
- reset()#
Resets a KSP context.
Collective.
Resets a KSP context to the
kspsetupcalled = 0
state and removes any allocated Vecs and Mats.See also
Source code at petsc4py/PETSc/KSP.pyx:1607
- Return type:
- setAppCtx(appctx)#
Set the optional user-defined context for the linear solver.
Not collective.
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
- 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.
Notes
Currently this option is not valid for all iterative methods.
- setComputeOperators(operators, args=None, kargs=None)#
Set routine to compute the linear operators.
Logically collective.
- Parameters:
- Return type:
Notes
The user provided function
Operators
will be called automatically at the very next call tosolve
. It will NOT be called at futuresolve
calls unless eithersetComputeOperators
orsetOperators
is called before thatsolve
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 eachsolve
.To reuse the same preconditioner for the next
solve
and not compute a new one based on the most recently computed matrix callKSPSetReusePreconditioner
.
- setComputeRHS(rhs, args=None, kargs=None)#
Set routine to compute the right-hand side of the linear system.
Logically collective.
- Parameters:
- Return type:
Notes
The routine you provide will be called each time you call
solve
to prepare the new right-hand side for that solve.See also
- 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.
Notes
Currently this option is not valid for all iterative methods.
- setConvergedReason(reason)#
Use
reason
property.Source code at petsc4py/PETSc/KSP.pyx:1871
- Parameters:
reason (ConvergedReason)
- Return type:
- 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:
- Return type:
Notes
If
length
is not provided orNone
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.See also
- setConvergenceTest(converged, args=None, kargs=None)#
Set the function to be used to determine convergence.
Logically collective.
- Parameters:
- Return type:
Notes
Must be called after the KSP type has been set so put this after a call to
setType
, orsetFromOptions
.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.
- setDM(dm)#
Set the
DM
that may be used by some preconditioners.Logically collective.
Notes
If this is used then the
KSP
will attempt to use theDM
to create the matrix and use the routine set withDM.setKSPComputeOperators
. UsesetDMActive(False)
to instead use the matrix you have provided withsetOperators
.A
DM
can only be used for solving one problem at a time because information about the problem is stored on theDM
, even when not using interfaces likeDM.setKSPComputeOperators
. UseDM.clone
to get a distinctDM
when solving different problems using the same function space.See also
PETSc.KSP
,DM
,DM.setKSPComputeOperators
,setOperators
,DM.clone
,KSPSetDM
- setDMActive(flag)#
DM
should be used to generate system matrix & RHS vector.Logically collective.
Notes
By default
setDM
sets theDM
as active, callsetDMActive(False)
aftersetDM(dm)
to not have theKSP
object use theDM
to generate the matrices.See also
- setErrorIfNotConverged(flag)#
Cause
solve
to generate an error if not converged.Logically collective.
See also
- 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.See also
Source code at petsc4py/PETSc/KSP.pyx:600
- Return type:
- setGMRESRestart(restart)#
Set number of iterations at which KSP restarts.
Logically collective.
Suitable KSPs are: KSPGMRES, KSPFGMRES and KSPLGMRES.
See also
- setHPDDMType(hpddm_type)#
Set the Krylov solver type.
Collective.
See also
- setInitialGuessKnoll(flag)#
Tell solver to use
PC.apply
to compute the initial guess.Logically collective.
This is the Knoll trick.
See also
- 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:
See also
- 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:
- Return type:
Notes
The default is to do nothing. To print the residual, or preconditioned residual if
setNormType(NORM_PRECONDITIONED)
was called, usemonitor
as the monitoring routine, with aPETSc.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.
- setNormType(normtype)#
Set the norm that is used for convergence testing.
Logically collective.
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.
- 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:
- Return type:
Notes
If you know the operator
A
has a null space you can useMat.setNullSpace
andMat.setTransposeNullSpace
to supply the null space toA
and theKSP
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
forA
orP
removes the matrix that is currently used.See also
- setOptionsPrefix(prefix)#
Set the prefix used for all
KSP
options in the database.Logically collective.
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 `
- setPC(pc)#
Set the preconditioner.
Collective.
Set the preconditioner to be used to calculate the application of the preconditioner on a vector.
- setPCSide(side)#
Set the preconditioning side.
Logically collective.
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.See also
PC.Side
, Working with PETSc options,getPCSide
,setNormType
,getNormType
,KSPSetPCSide
- setPostSolve(postsolve, args=None, kargs=None)#
Set the function that is called at the end of each
KSP.solve
.Logically collective.
- Parameters:
- Return type:
See also
- setPreSolve(presolve, args=None, kargs=None)#
Set the function that is called at the beginning of each
KSP.solve
.Logically collective.
- Parameters:
- Return type:
See also
- setPythonContext(context=None)#
Set the instance of the class implementing Python methods.
Not collective.
- setPythonType(py_type)#
Set the fully qualified Python name of the class to be used.
Collective.
- 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. Or
DETERMINE
to use the value when the object’s type was set.atol (float | None) – The absolute convergence tolerance absolute size of the (possibly preconditioned) residual norm. Or
DETERMINE
to use the value when the object’s type was set.dtol – The divergence tolerance, amount (possibly preconditioned) residual norm can increase before
KSPConvergedDefault
concludes that the method is diverging. OrDETERMINE
to use the value when the object’s type was set.max_it (int | None) – Maximum number of iterations to use. Or
DETERMINE
to use the value when the object’s type was set.
- Return type:
Notes
Use
None
to retain the default value of any of the tolerances.
- setType(ksp_type)#
Build the
KSP
data structure for a particularType
.Logically collective.
Notes
See
Type
for available methods (for instance,Type.CG
orType.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
- setUp()#
Set up internal data structures for an iterative solver.
Collective.
See also
Source code at petsc4py/PETSc/KSP.pyx:1595
- Return type:
- 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
Source code at petsc4py/PETSc/KSP.pyx:1622
- Return type:
- 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:
- Return type:
See also
- solve(b, x)#
Solve the linear system.
Collective.
Notes
If one uses
setDM
thenx
orb
need not be passed. UsegetSolution
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. CallgetConvergedReason
to determine if the solver converged or failed and why. The option-ksp_error_if_not_converged
or functionsetErrorIfNotConverged
will causesolve
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
andMat.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 asPC.Type.LU
orPC.Type.ILU
, then its=1. SeesetTolerances
for more details.Understanding Convergence
The routines
setMonitor
andcomputeEigenvalues
provide information on additional options to monitor convergence and print eigenvalue information.
- solveTranspose(b, x)#
Solve the transpose of a linear system.
Collective.
Notes
For complex numbers this solve the non-Hermitian transpose system.
See also
- view(viewer=None)#
Print the KSP data structure.
Collective.
See also
Attributes Documentation
- appctx#
The solver application context.
- atol#
The absolute tolerance of the solver.
- divtol#
The divergence tolerance of the solver.
- 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_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.