xref: /petsc/src/binding/petsc4py/src/petsc4py/PETSc/KSP.pyx (revision 3220ff8572602716d60bdda8b51773ebceb3c8ea)
1# --------------------------------------------------------------------
2
3class KSPType(object):
4    """KSP Type.
5
6    The available types are:
7
8    `RICHARDSON`
9        The preconditioned Richardson iterative method
10        `petsc.KSPRICHARDSON`.
11    `CHEBYSHEV`
12        The preconditioned Chebyshev iterative method.
13        `petsc.KSPCHEBYSHEV`.
14    `CG`
15        The Preconditioned Conjugate Gradient (PCG) iterative method.
16        `petsc.KSPCG`
17    `GROPPCG`
18        A pipelined conjugate gradient method (Gropp).
19        `petsc.KSPGROPPCG`
20    `PIPECG`
21        A pipelined conjugate gradient method.
22        `petsc.KSPPIPECG`
23    `PIPECGRR`
24        Pipelined Conjugate Gradients with Residual Replacement.
25        `petsc.KSPPIPECGRR`
26    `PIPELCG`
27        Deep pipelined (length l) Conjugate Gradient method.
28        `petsc.KSPPIPELCG`
29    `PIPEPRCG`
30        Pipelined predict-and-recompute conjugate gradient method.
31        `petsc.KSPPIPEPRCG`
32    `PIPECG2`
33        Pipelined conjugate gradient method with a single non-blocking
34        reduction per two iterations. `petsc.KSPPIPECG2`
35    `CGNE`
36        Applies the preconditioned conjugate gradient method to the
37        normal equations without explicitly forming AᵀA. `petsc.KSPCGNE`
38    `NASH`
39        Conjugate gradient method subject to a constraint
40        on the solution norm. `petsc.KSPNASH`
41    `STCG`
42        Conjugate gradient method subject to a constraint on the
43        solution norm. `petsc.KSPSTCG`
44    `GLTR`
45        Conjugate gradient method subject to a constraint on the
46        solution norm. `petsc.KSPGLTR`
47    `FCG`
48        Flexible Conjugate Gradient method (FCG). Unlike most KSP
49        methods this allows the preconditioner to be nonlinear.
50        `petsc.KSPFCG`
51    `PIPEFCG`
52        Pipelined, Flexible Conjugate Gradient method.
53        `petsc.KSPPIPEFCG`
54    `GMRES`
55        Generalized Minimal Residual method with restart.
56        `petsc.KSPGMRES`
57    `PIPEFGMRES`
58        Pipelined (1-stage) Flexible Generalized Minimal Residual
59        method. `petsc.KSPPIPEFGMRES`
60    `FGMRES`
61        Implements the Flexible Generalized Minimal Residual method.
62        `petsc.KSPFGMRES`
63    `LGMRES`
64        Augments the standard Generalized Minimal Residual method
65        approximation space with approximations to the error from
66        previous restart cycles. `petsc.KSPLGMRES`
67    `DGMRES`
68        Deflated Generalized Minimal Residual method. In this
69        implementation, the adaptive strategy allows to switch to the
70        deflated GMRES when the stagnation occurs. `petsc.KSPDGMRES`
71    `PGMRES`
72        Pipelined Generalized Minimal Residual method.
73        `petsc.KSPPGMRES`
74    `TCQMR`
75        A variant of Quasi Minimal Residual (QMR).
76        `petsc.KSPTCQMR`
77    `BCGS`
78        Stabilized version of Biconjugate Gradient (BiCGStab) method.
79        `petsc.KSPBCGS`
80    `IBCGS`
81        Improved Stabilized version of BiConjugate Gradient (IBiCGStab)
82        method in an alternative form to have only a single global
83        reduction operation instead of the usual 3 (or 4).
84        `petsc.KSPIBCGS`
85    `QMRCGS`
86        Quasi- Minimal Residual variant of the Bi-CGStab algorithm
87        (QMRCGStab) method. `petsc.KSPQMRCGS`
88    `FBCGS`
89        Flexible Stabilized version of BiConjugate Gradient (BiCGStab)
90        method. `petsc.KSPFBCGS`
91    `FBCGSR`
92        A mathematically equivalent variant of flexible stabilized
93        BiConjugate Gradient (BiCGStab). `petsc.KSPFBCGSR`
94    `BCGSL`
95        Variant of the L-step stabilized BiConjugate Gradient
96        (BiCGStab(L)) algorithm. Uses "L-step" Minimal Residual (MR)
97        polynomials. The variation concerns cases when some parameters
98        are negative due to round-off. `petsc.KSPBCGSL`
99    `PIPEBCGS`
100        Pipelined stabilized BiConjugate Gradient (BiCGStab) method.
101        `petsc.KSPPIPEBCGS`
102    `CGS`
103        Conjugate Gradient Squared method.
104        `petsc.KSPCGS`
105    `TFQMR`
106        A Transpose Tree Quasi- Minimal Residual (QMR).
107        `petsc.KSPCR`
108    `CR`
109        (Preconditioned) Conjugate Residuals (CR) method.
110        `petsc.KSPCR`
111    `PIPECR`
112        Pipelined Conjugate Residual (CR) method.
113        `petsc.KSPPIPECR`
114    `LSQR`
115        Least squares solver.
116        `petsc.KSPLSQR`
117    `PREONLY`
118        Applies ONLY the preconditioner exactly once. This may be used
119        in inner iterations, where it is desired to allow multiple
120        iterations as well as the "0-iteration" case. It is commonly
121        used with the direct solver preconditioners like PCLU and
122        PCCHOLESKY. There is an alias of KSPNONE.
123        `petsc.KSPPREONLY`
124    `NONE`
125        No solver
126        ``KSPNONE``
127    `QCG`
128        Conjugate Gradient (CG) method subject to a constraint on the
129        solution norm. `petsc.KSPQCG`
130    `BICG`
131        Implements the Biconjugate gradient method (BiCG).
132        Similar to running the conjugate gradient on the normal equations.
133        `petsc.KSPBICG`
134    `MINRES`
135        Minimum Residual (MINRES) method.
136        `petsc.KSPMINRES`
137    `SYMMLQ`
138        Symmetric LQ method (SymmLQ). Uses LQ decomposition (lower
139        trapezoidal).
140        `petsc.KSPSYMMLQ`
141    `LCD`
142        Left Conjugate Direction (LCD) method.
143        `petsc.KSPLCD`
144    `PYTHON`
145        Python shell solver. Call Python function to implement solver.
146        ``KSPPYTHON``
147    `GCR`
148        Preconditioned flexible Generalized Conjugate Residual (GCR)
149        method.
150        `petsc.KSPGCR`
151    `PIPEGCR`
152        Pipelined Generalized Conjugate Residual method.
153        `petsc.KSPPIPEGCR`
154    `TSIRM`
155        Two-Stage Iteration with least-squares Residual Minimization
156        method. `petsc.KSPTSIRM`
157    `CGLS`
158        Conjugate Gradient method for Least-Squares problems. Supports
159        non-square (rectangular) matrices. `petsc.KSPCGLS`
160    `FETIDP`
161        Dual-Primal (DP) Finite Element Tearing and Interconnect (FETI)
162        method. `petsc.KSPFETIDP`
163    `HPDDM`
164        Interface with the HPDDM library. This KSP may be used to
165        further select methods that are currently not implemented
166        natively in PETSc, e.g., GCRODR, a recycled Krylov
167        method which is similar to KSPLGMRES. `petsc.KSPHPDDM`
168
169    See Also
170    --------
171    petsc_options, petsc.KSPType
172
173    """
174    RICHARDSON = S_(KSPRICHARDSON)
175    CHEBYSHEV  = S_(KSPCHEBYSHEV)
176    CG         = S_(KSPCG)
177    GROPPCG    = S_(KSPGROPPCG)
178    PIPECG     = S_(KSPPIPECG)
179    PIPECGRR   = S_(KSPPIPECGRR)
180    PIPELCG    = S_(KSPPIPELCG)
181    PIPEPRCG   = S_(KSPPIPEPRCG)
182    PIPECG2    = S_(KSPPIPECG2)
183    CGNE       = S_(KSPCGNE)
184    NASH       = S_(KSPNASH)
185    STCG       = S_(KSPSTCG)
186    GLTR       = S_(KSPGLTR)
187    FCG        = S_(KSPFCG)
188    PIPEFCG    = S_(KSPPIPEFCG)
189    GMRES      = S_(KSPGMRES)
190    PIPEFGMRES = S_(KSPPIPEFGMRES)
191    FGMRES     = S_(KSPFGMRES)
192    LGMRES     = S_(KSPLGMRES)
193    DGMRES     = S_(KSPDGMRES)
194    PGMRES     = S_(KSPPGMRES)
195    TCQMR      = S_(KSPTCQMR)
196    BCGS       = S_(KSPBCGS)
197    IBCGS      = S_(KSPIBCGS)
198    QMRCGS     = S_(KSPQMRCGS)
199    FBCGS      = S_(KSPFBCGS)
200    FBCGSR     = S_(KSPFBCGSR)
201    BCGSL      = S_(KSPBCGSL)
202    PIPEBCGS   = S_(KSPPIPEBCGS)
203    CGS        = S_(KSPCGS)
204    TFQMR      = S_(KSPTFQMR)
205    CR         = S_(KSPCR)
206    PIPECR     = S_(KSPPIPECR)
207    LSQR       = S_(KSPLSQR)
208    PREONLY    = S_(KSPPREONLY)
209    NONE       = S_(KSPNONE)
210    QCG        = S_(KSPQCG)
211    BICG       = S_(KSPBICG)
212    MINRES     = S_(KSPMINRES)
213    SYMMLQ     = S_(KSPSYMMLQ)
214    LCD        = S_(KSPLCD)
215    PYTHON     = S_(KSPPYTHON)
216    GCR        = S_(KSPGCR)
217    PIPEGCR    = S_(KSPPIPEGCR)
218    TSIRM      = S_(KSPTSIRM)
219    CGLS       = S_(KSPCGLS)
220    FETIDP     = S_(KSPFETIDP)
221    HPDDM      = S_(KSPHPDDM)
222
223
224class KSPNormType(object):
225    """KSP norm type.
226
227    The available norm types are:
228
229    `NONE`
230        Skips computing the norm, this should generally only be used if
231        you are using the Krylov method as a smoother with a fixed
232        small number of iterations. Implicitly sets
233        `petsc.KSPConvergedSkip` as KSP convergence test. Note that
234        certain algorithms such as `Type.GMRES` ALWAYS require the norm
235        calculation, for these methods the norms are still computed,
236        they are just not used in the convergence test.
237    `PRECONDITIONED`
238        The default for left preconditioned solves, uses the l₂ norm of
239        the preconditioned residual P⁻¹(b - Ax).
240    `UNPRECONDITIONED`
241        Uses the l₂ norm of the true b - Ax residual.
242    `NATURAL`
243        Supported by `Type.CG`, `Type.CR`, `Type.CGNE`, `Type.CGS`.
244
245    """
246    # native
247    NORM_DEFAULT          = KSP_NORM_DEFAULT
248    NORM_NONE             = KSP_NORM_NONE
249    NORM_PRECONDITIONED   = KSP_NORM_PRECONDITIONED
250    NORM_UNPRECONDITIONED = KSP_NORM_UNPRECONDITIONED
251    NORM_NATURAL          = KSP_NORM_NATURAL
252    # aliases
253    DEFAULT          = NORM_DEFAULT
254    NONE = NO        = NORM_NONE
255    PRECONDITIONED   = NORM_PRECONDITIONED
256    UNPRECONDITIONED = NORM_UNPRECONDITIONED
257    NATURAL          = NORM_NATURAL
258
259
260class KSPConvergedReason(object):
261    """KSP Converged Reason.
262
263    `CONVERGED_ITERATING`
264        Still iterating
265    `ITERATING`
266        Still iterating
267
268    `CONVERGED_RTOL_NORMAL_EQUATIONS`
269        Undocumented.
270    `CONVERGED_ATOL_NORMAL_EQUATIONS`
271        Undocumented.
272    `CONVERGED_RTOL`
273        ∥r∥ <= rtolnorm(b) or rtolnorm(b - Ax₀)
274    `CONVERGED_ATOL`
275        ∥r∥ <= atol
276    `CONVERGED_ITS`
277        Used by the `Type.PREONLY` solver after the single iteration of the
278        preconditioner is applied. Also used when the
279        `petsc.KSPConvergedSkip` convergence test routine is set in KSP.
280    `CONVERGED_NEG_CURVE`
281        Undocumented.
282    `CONVERGED_STEP_LENGTH`
283        Undocumented.
284    `CONVERGED_HAPPY_BREAKDOWN`
285        Undocumented.
286
287    `DIVERGED_NULL`
288        Undocumented.
289    `DIVERGED_MAX_IT`
290        Ran out of iterations before any convergence criteria was
291        reached.
292    `DIVERGED_DTOL`
293        norm(r) >= dtol*norm(b)
294    `DIVERGED_BREAKDOWN`
295        A breakdown in the Krylov method was detected so the method
296        could not continue to enlarge the Krylov space. Could be due to
297        a singular matrix or preconditioner. In KSPHPDDM, this is also
298        returned when some search directions within a block are
299        collinear.
300    `DIVERGED_BREAKDOWN_BICG`
301        A breakdown in the KSPBICG method was detected so the method
302        could not continue to enlarge the Krylov space.
303    `DIVERGED_NONSYMMETRIC`
304        It appears the operator or preconditioner is not symmetric and
305        this Krylov method (`Type.CG`, `Type.MINRES`, `Type.CR`)
306        requires symmetry.
307    `DIVERGED_INDEFINITE_PC`
308        It appears the preconditioner is indefinite (has both positive
309        and negative eigenvalues) and this Krylov method (`Type.CG`)
310        requires it to be positive definite.
311    `DIVERGED_NANORINF`
312        Undocumented.
313    `DIVERGED_INDEFINITE_MAT`
314        Undocumented.
315    `DIVERGED_PCSETUP_FAILED`
316        It was not possible to build or use the requested
317        preconditioner. This is usually due to a zero pivot in a
318        factorization. It can also result from a failure in a
319        subpreconditioner inside a nested preconditioner such as
320        `PC.Type.FIELDSPLIT`.
321
322    See Also
323    --------
324    `petsc.KSPConvergedReason`
325
326    """
327    # iterating
328    CONVERGED_ITERATING       = KSP_CONVERGED_ITERATING
329    ITERATING                 = KSP_CONVERGED_ITERATING
330    # converged
331    CONVERGED_RTOL_NORMAL_EQUATIONS = KSP_CONVERGED_RTOL_NORMAL_EQUATIONS
332    CONVERGED_ATOL_NORMAL_EQUATIONS = KSP_CONVERGED_ATOL_NORMAL_EQUATIONS
333    CONVERGED_RTOL            = KSP_CONVERGED_RTOL
334    CONVERGED_ATOL            = KSP_CONVERGED_ATOL
335    CONVERGED_ITS             = KSP_CONVERGED_ITS
336    CONVERGED_NEG_CURVE       = KSP_CONVERGED_NEG_CURVE
337    CONVERGED_STEP_LENGTH     = KSP_CONVERGED_STEP_LENGTH
338    CONVERGED_HAPPY_BREAKDOWN = KSP_CONVERGED_HAPPY_BREAKDOWN
339    # diverged
340    DIVERGED_NULL             = KSP_DIVERGED_NULL
341    DIVERGED_MAX_IT           = KSP_DIVERGED_MAX_IT
342    DIVERGED_DTOL             = KSP_DIVERGED_DTOL
343    DIVERGED_BREAKDOWN        = KSP_DIVERGED_BREAKDOWN
344    DIVERGED_BREAKDOWN_BICG   = KSP_DIVERGED_BREAKDOWN_BICG
345    DIVERGED_NONSYMMETRIC     = KSP_DIVERGED_NONSYMMETRIC
346    DIVERGED_INDEFINITE_PC    = KSP_DIVERGED_INDEFINITE_PC
347    DIVERGED_NANORINF         = KSP_DIVERGED_NANORINF
348    DIVERGED_INDEFINITE_MAT   = KSP_DIVERGED_INDEFINITE_MAT
349    DIVERGED_PCSETUP_FAILED   = KSP_DIVERGED_PC_FAILED
350
351
352class KSPHPDDMType(object):
353    """The *HPDDM* Krylov solver type."""
354    GMRES                     = KSP_HPDDM_TYPE_GMRES
355    BGMRES                    = KSP_HPDDM_TYPE_BGMRES
356    CG                        = KSP_HPDDM_TYPE_CG
357    BCG                       = KSP_HPDDM_TYPE_BCG
358    GCRODR                    = KSP_HPDDM_TYPE_GCRODR
359    BGCRODR                   = KSP_HPDDM_TYPE_BGCRODR
360    BFBCG                     = KSP_HPDDM_TYPE_BFBCG
361    PREONLY                   = KSP_HPDDM_TYPE_PREONLY
362
363
364class KSPDMActive(object):
365    """Argument to KSPSetDMActive()."""
366    OPERATOR      = KSP_DMACTIVE_OPERATOR
367    RHS           = KSP_DMACTIVE_RHS
368    INITIAL_GUESS = KSP_DMACTIVE_INITIAL_GUESS
369    ALL           = KSP_DMACTIVE_ALL
370
371
372# --------------------------------------------------------------------
373
374
375cdef class KSP(Object):
376    """Abstract PETSc object that manages all Krylov methods.
377
378    This is the object that manages the linear solves in PETSc (even
379    those such as direct solvers that do no use Krylov accelerators).
380
381    Notes
382    -----
383    When a direct solver is used, but no Krylov solver is used, the KSP
384    object is still used but with a `Type.PREONLY`, meaning that
385    only application of the preconditioner is used as the linear
386    solver.
387
388    See Also
389    --------
390    create, setType, SNES, TS, PC, Type.CG, Type.GMRES,
391    petsc.KSP
392
393    """
394
395    Type            = KSPType
396    NormType        = KSPNormType
397    ConvergedReason = KSPConvergedReason
398    HPDDMType       = KSPHPDDMType
399    DMActive        = KSPDMActive
400
401    # --- xxx ---
402
403    def __cinit__(self):
404        self.obj = <PetscObject*> &self.ksp
405        self.ksp = NULL
406
407    def __call__(self, Vec b, Vec x = None) -> Vec:
408        """Solve linear system.
409
410        Collective.
411
412        Parameters
413        ----------
414        b
415            Right hand side vector.
416        x
417            Solution vector.
418
419        Notes
420        -----
421        Shortcut for `solve`, which returns the solution vector.
422
423        See Also
424        --------
425        solve, petsc_options, petsc.KSPSolve
426
427        """
428        if x is None: # XXX do this better
429            x = self.getOperators()[0].createVecLeft()
430        self.solve(b, x)
431        return x
432
433    # --- xxx ---
434
435    def view(self, Viewer viewer=None) -> None:
436        """Print the KSP data structure.
437
438        Collective.
439
440        Parameters
441        ----------
442        viewer
443            Viewer used to display the KSP.
444
445        See Also
446        --------
447        petsc.KSPView
448
449        """
450        cdef PetscViewer vwr = NULL
451        if viewer is not None: vwr = viewer.vwr
452        CHKERR(KSPView(self.ksp, vwr))
453
454    def destroy(self) -> Self:
455        """Destroy KSP context.
456
457        Collective.
458
459        See Also
460        --------
461        petsc.KSPDestroy
462
463        """
464        CHKERR(KSPDestroy(&self.ksp))
465        return self
466
467    def create(self, comm: Comm | None = None) -> Self:
468        """Create the KSP context.
469
470        Collective.
471
472        See Also
473        --------
474        petsc.KSPCreate
475
476        """
477        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_DEFAULT)
478        cdef PetscKSP newksp = NULL
479        CHKERR(KSPCreate(ccomm, &newksp))
480        CHKERR(PetscCLEAR(self.obj)); self.ksp = newksp
481        return self
482
483    def setType(self, ksp_type: Type | str) -> None:
484        """Build the `KSP` data structure for a particular `Type`.
485
486        Logically collective.
487
488        Parameters
489        ----------
490        ksp_type
491            KSP Type object
492
493        Notes
494        -----
495        See `Type` for available methods (for instance, `Type.CG` or
496        `Type.GMRES`).
497
498        Normally, it is best to use the `setFromOptions` command
499        and then set the KSP type from the options database rather than
500        by using this routine. Using the options database provides the
501        user with maximum flexibility in evaluating the many different
502        Krylov methods. This method is provided for those situations
503        where it is necessary to set the iterative solver independently
504        of the command line or options database. This might be the
505        case, for example, when the choice of iterative solver changes
506        during the execution of the program, and the user's application
507        is taking responsibility for choosing the appropriate method.
508        In other words, this routine is not for beginners.
509
510        See Also
511        --------
512        petsc.KSPSetType
513
514        """
515        cdef PetscKSPType cval = NULL
516        ksp_type = str2bytes(ksp_type, &cval)
517        CHKERR(KSPSetType(self.ksp, cval))
518
519    def getType(self) -> str:
520        """Return the KSP type as a string from the `KSP` object.
521
522        Not collective.
523
524        See Also
525        --------
526        petsc.KSPGetType
527
528        """
529        cdef PetscKSPType cval = NULL
530        CHKERR(KSPGetType(self.ksp, &cval))
531        return bytes2str(cval)
532
533    def setOptionsPrefix(self, prefix: str | None) -> None:
534        """Set the prefix used for all `KSP` options in the database.
535
536        Logically collective.
537
538        Parameters
539        ----------
540        prefix
541            The options prefix.
542
543        Notes
544        -----
545        A hyphen (-) must NOT be given at the beginning of the prefix
546        name. The first character of all runtime options is
547        AUTOMATICALLY the hyphen. For example, to distinguish between
548        the runtime options for two different `KSP` contexts, one could
549        call
550        ```
551        KSPSetOptionsPrefix(ksp1, "sys1_")
552        KSPSetOptionsPrefix(ksp2, "sys2_")
553        ```
554
555        This would enable use of different options for each system,
556        such as
557        ```
558        -sys1_ksp_type gmres -sys1_ksp_rtol 1.e-3
559        -sys2_ksp_type bcgs  -sys2_ksp_rtol 1.e-4
560        ```
561
562        See Also
563        --------
564        petsc_options, petsc.KSPSetOptionsPrefix
565
566        """
567        cdef const char *cval = NULL
568        prefix = str2bytes(prefix, &cval)
569        CHKERR(KSPSetOptionsPrefix(self.ksp, cval))
570
571    def getOptionsPrefix(self) -> str:
572        """Return the prefix used for all `KSP` options in the database.
573
574        Not collective.
575
576        See Also
577        --------
578        petsc.KSPGetOptionsPrefix
579
580        """
581        cdef const char *cval = NULL
582        CHKERR(KSPGetOptionsPrefix(self.ksp, &cval))
583        return bytes2str(cval)
584
585    def appendOptionsPrefix(self, prefix: str | None) -> None:
586        """Append to prefix used for all `KSP` options in the database.
587
588        Logically collective.
589
590        Parameters
591        ----------
592        prefix
593            The options prefix to append.
594
595        Notes
596        -----
597        A hyphen (-) must NOT be given at the beginning of the prefix
598        name. The first character of all runtime options is
599        AUTOMATICALLY the hyphen.
600
601        See Also
602        --------
603        petsc.KSPAppendOptionsPrefix
604
605        """
606        cdef const char *cval = NULL
607        prefix = str2bytes(prefix, &cval)
608        CHKERR(KSPAppendOptionsPrefix(self.ksp, cval))
609
610    def setFromOptions(self) -> None:
611        """Set `KSP` options from the options database.
612
613        Collective.
614
615        This routine must be called before `setUp` if the user is
616        to be allowed to set the Krylov type.
617
618        See Also
619        --------
620        petsc_options, petsc.KSPSetFromOptions
621
622        """
623        CHKERR(KSPSetFromOptions(self.ksp))
624
625    # --- application context ---
626
627    def setAppCtx(self, appctx: Any) -> None:
628        """Set the optional user-defined context for the linear solver.
629
630        Not collective.
631
632        Parameters
633        ----------
634        appctx
635            The user defined context
636
637        Notes
638        -----
639        The user context is a way for users to attach any information
640        to the `KSP` that they may need later when interacting with
641        the solver.
642
643        See Also
644        --------
645        getAppCtx
646
647        """
648        self.set_attr('__appctx__', appctx)
649
650    def getAppCtx(self) -> Any:
651        """Return the user-defined context for the linear solver.
652
653        Not collective.
654
655        See Also
656        --------
657        setAppCtx
658
659        """
660        return self.get_attr('__appctx__')
661
662    # --- discretization space ---
663
664    def getDM(self) -> DM:
665        """Return the `DM` that may be used by some preconditioners.
666
667        Not collective.
668
669        See Also
670        --------
671        PETSc.KSP, DM, petsc.KSPGetDM
672
673        """
674        cdef PetscDM newdm = NULL
675        CHKERR(KSPGetDM(self.ksp, &newdm))
676        cdef DM dm = subtype_DM(newdm)()
677        dm.dm = newdm
678        CHKERR(PetscINCREF(dm.obj))
679        return dm
680
681    def setDM(self, DM dm) -> None:
682        """Set the `DM` that may be used by some preconditioners.
683
684        Logically collective.
685
686        Parameters
687        ----------
688        dm
689            The `DM` object, cannot be `None`.
690
691        Notes
692        -----
693        If this is used then the `KSP` will attempt to use the `DM` to
694        create the matrix and use the routine set with
695        `DM.setKSPComputeOperators`. Use
696        ``setDMActive(KSP.DMActive.OPERATOR, False)``
697        to instead use the matrix you have provided with
698        `setOperators`.
699
700        A `DM` can only be used for solving one problem at a time
701        because information about the problem is stored on the `DM`,
702        even when not using interfaces like
703        `DM.setKSPComputeOperators`. Use `DM.clone` to get a distinct
704        `DM` when solving different problems using the same function
705        space.
706
707        See Also
708        --------
709        PETSc.KSP, DM, DM.setKSPComputeOperators, setOperators, DM.clone
710        petsc.KSPSetDM
711
712        """
713        CHKERR(KSPSetDM(self.ksp, dm.dm))
714
715    def setDMActive(self, dmactive: DMActive, flag: bool) -> None:
716        """`DM` should be used to generate system matrix & RHS vector.
717
718        Logically collective.
719
720        Parameters
721        ----------
722        flag
723            Whether to use the `DM`.
724
725        Notes
726        -----
727        By default `setDM` sets the `DM` as active, call
728        ``setDMActive(KSP.DMactive.ALL, False)`` after ``setDM(dm)`` to not
729        have the `KSP` object use the `DM` to generate the matrices.
730
731        See Also
732        --------
733        PETSc.KSP, DM, setDM, petsc.KSPSetDMActive
734
735        """
736        cdef PetscBool cflag = asBool(flag)
737        CHKERR(KSPSetDMActive(self.ksp, dmactive, cflag))
738
739    # --- operators and preconditioner ---
740
741    def setComputeRHS(
742        self,
743        rhs: KSPRHSFunction,
744        args: tuple[Any, ...] | None = None,
745        kargs: dict[str, Any] | None = None) -> None:
746        """Set routine to compute the right-hand side of the linear system.
747
748        Logically collective.
749
750        Parameters
751        ----------
752        rhs
753            Function which computes the right-hand side.
754        args
755            Positional arguments for callback function ``rhs``.
756        kargs
757            Keyword arguments for callback function ``rhs``.
758
759        Notes
760        -----
761        The routine you provide will be called each time you call `solve`
762        to prepare the new right-hand side for that solve.
763
764        See Also
765        --------
766        PETSc.KSP, solve, petsc.KSPSetComputeRHS
767
768        """
769        if args  is None: args  = ()
770        if kargs is None: kargs = {}
771        context = (rhs, args, kargs)
772        self.set_attr('__rhs__', context)
773        CHKERR(KSPSetComputeRHS(self.ksp, KSP_ComputeRHS, <void*>context))
774
775    def setComputeOperators(
776        self,
777        operators: KSPOperatorsFunction,
778        args: tuple[Any, ...] | None = None,
779        kargs: dict[str, Any] | None = None) -> None:
780        """Set routine to compute the linear operators.
781
782        Logically collective.
783
784        Parameters
785        ----------
786        operators
787            Function which computes the operators.
788        args
789            Positional arguments for callback function ``operators``.
790        kargs
791            Keyword arguments for callback function ``operators``.
792
793        Notes
794        -----
795        The user provided function `operators` will be called
796        automatically at the very next call to `solve`. It will NOT
797        be called at future `solve` calls unless either
798        `setComputeOperators` or `setOperators` is called
799        before that `solve` is called. This allows the same system
800        to be solved several times with different right-hand side
801        functions, but is a confusing API since one might expect it to
802        be called for each `solve`.
803
804        To reuse the same preconditioner for the next `solve` and
805        not compute a new one based on the most recently computed
806        matrix call `petsc.KSPSetReusePreconditioner`.
807
808        See Also
809        --------
810        PETSc.KSP, solve, setOperators, petsc.KSPSetComputeOperators
811        petsc.KSPSetReusePreconditioner
812
813        """
814        if args  is None: args  = ()
815        if kargs is None: kargs = {}
816        context = (operators, args, kargs)
817        self.set_attr('__operators__', context)
818        CHKERR(KSPSetComputeOperators(self.ksp, KSP_ComputeOps, <void*>context))
819
820    def setOperators(self, Mat A=None, Mat P=None) -> None:
821        """Set matrix associated with the linear system.
822
823        Collective.
824
825        Set the matrix associated with the linear system and a
826        (possibly) different one from which the preconditioner will be
827        built.
828
829        Parameters
830        ----------
831        A
832            Matrix that defines the linear system.
833        P
834            Matrix to be used in constructing the preconditioner,
835            usually the same as ``A``.
836
837        Notes
838        -----
839        This is equivalent to ``pc = ksp.getPC(); pc.setOperators(A, P)``
840        but is the preferred approach.
841
842        If you know the operator ``A`` has a null space you can use
843        `Mat.setNullSpace` and `Mat.setTransposeNullSpace` to supply the
844        null space to ``A`` and the `KSP` solvers will automatically use
845        that null space as needed during the solution process.
846
847        All future calls to `setOperators` must use the same size
848        matrices!
849
850        Passing `None` for ``A`` or ``P`` removes the matrix that is
851        currently used.
852
853        See Also
854        --------
855        PETSc.KSP, solve, setComputeOperators, petsc.KSPSetOperators
856
857        """
858        cdef PetscMat amat=NULL
859        if A is not None: amat = A.mat
860        cdef PetscMat pmat=amat
861        if P is not None: pmat = P.mat
862        CHKERR(KSPSetOperators(self.ksp, amat, pmat))
863
864    def getOperators(self) -> tuple[Mat, Mat]:
865        """Return the matrix associated with the linear system.
866
867        Collective.
868
869        Return the matrix associated with the linear system and a
870        (possibly) different one used to construct the preconditioner.
871
872        Returns
873        -------
874        A : Mat
875            Matrix that defines the linear system.
876        P : Mat
877            Matrix to be used in constructing the preconditioner,
878            usually the same as ``A``.
879
880        See Also
881        --------
882        PETSc.KSP, solve, setOperators, petsc.KSPGetOperators
883
884        """
885        cdef Mat A = Mat(), P = Mat()
886        CHKERR(KSPGetOperators(self.ksp, &A.mat, &P.mat))
887        CHKERR(PetscINCREF(A.obj))
888        CHKERR(PetscINCREF(P.obj))
889        return (A, P)
890
891    def setPC(self, PC pc) -> None:
892        """Set the preconditioner.
893
894        Collective.
895
896        Set the preconditioner to be used to calculate the application
897        of the preconditioner on a vector.
898
899        Parameters
900        ----------
901        pc
902            The preconditioner object
903
904        See Also
905        --------
906        PETSc.KSP, getPC, petsc.KSPSetPC
907
908        """
909        CHKERR(KSPSetPC(self.ksp, pc.pc))
910
911    def getPC(self) -> PC:
912        """Return the preconditioner.
913
914        Not collective.
915
916        See Also
917        --------
918        PETSc.KSP, setPC, petsc.KSPGetPC
919
920        """
921        cdef PC pc = PC()
922        CHKERR(KSPGetPC(self.ksp, &pc.pc))
923        CHKERR(PetscINCREF(pc.obj))
924        return pc
925
926    # --- tolerances and convergence ---
927
928    def setTolerances(
929        self,
930        rtol: float | None = None,
931        atol: float | None = None,
932        divtol: float | None = None,
933        max_it: int | None = None) -> None:
934        """Set various tolerances used by the KSP convergence testers.
935
936        Logically collective.
937
938        Set the relative, absolute, divergence, and maximum iteration
939        tolerances used by the default KSP convergence testers.
940
941        Parameters
942        ----------
943        rtol
944            The relative convergence tolerance, relative decrease in
945            the (possibly preconditioned) residual norm.
946            Or `DETERMINE` to use the value when
947            the object's type was set.
948        atol
949            The absolute convergence tolerance absolute size of the
950            (possibly preconditioned) residual norm.
951            Or `DETERMINE` to use the value when
952            the object's type was set.
953        dtol
954            The divergence tolerance, amount (possibly preconditioned)
955            residual norm can increase before
956            `petsc.KSPConvergedDefault` concludes that the method is
957            diverging.
958            Or `DETERMINE` to use the value when
959            the object's type was set.
960        max_it
961            Maximum number of iterations to use.
962            Or `DETERMINE` to use the value when
963            the object's type was set.
964
965        Notes
966        -----
967        Use `None` to retain the default value of any of the
968        tolerances.
969
970        See Also
971        --------
972        petsc_options, getTolerances, setConvergenceTest
973        petsc.KSPSetTolerances, petsc.KSPConvergedDefault
974
975        """
976        cdef PetscReal crtol, catol, cdivtol
977        crtol = catol = cdivtol = PETSC_CURRENT
978        if rtol   is not None: crtol   = asReal(rtol)
979        if atol   is not None: catol   = asReal(atol)
980        if divtol is not None: cdivtol = asReal(divtol)
981        cdef PetscInt cmaxits = PETSC_CURRENT
982        if max_it is not None: cmaxits = asInt(max_it)
983        CHKERR(KSPSetTolerances(self.ksp, crtol, catol, cdivtol, cmaxits))
984
985    def getTolerances(self) -> tuple[float, float, float, int]:
986        """Return various tolerances used by the KSP convergence tests.
987
988        Not collective.
989
990        Return the relative, absolute, divergence, and maximum iteration
991        tolerances used by the default KSP convergence tests.
992
993        Returns
994        -------
995        rtol : float
996            The relative convergence tolerance
997        atol : float
998            The absolute convergence tolerance
999        dtol : float
1000            The divergence tolerance
1001        maxits : int
1002            Maximum number of iterations
1003
1004        See Also
1005        --------
1006        setTolerances, petsc.KSPGetTolerances
1007
1008        """
1009        cdef PetscReal crtol=0, catol=0, cdivtol=0
1010        cdef PetscInt cmaxits=0
1011        CHKERR(KSPGetTolerances(self.ksp, &crtol, &catol, &cdivtol, &cmaxits))
1012        return (toReal(crtol), toReal(catol), toReal(cdivtol), toInt(cmaxits))
1013
1014    def setConvergenceTest(
1015        self,
1016        converged: KSPConvergenceTestFunction,
1017        args: tuple[Any, ...] | None = None,
1018        kargs: dict[str, Any] | None = None) -> None:
1019        """Set the function to be used to determine convergence.
1020
1021        Logically collective.
1022
1023        Parameters
1024        ----------
1025        converged
1026            Callback which computes the convergence.
1027        args
1028            Positional arguments for callback function.
1029        kargs
1030            Keyword arguments for callback function.
1031
1032        Notes
1033        -----
1034        Must be called after the KSP type has been set so put this
1035        after a call to `setType`, or `setFromOptions`.
1036
1037        The default is a combination of relative and absolute
1038        tolerances. The residual value that is tested may be an
1039        approximation; routines that need exact values should compute
1040        them.
1041
1042        See Also
1043        --------
1044        addConvergenceTest, ConvergedReason, setTolerances,
1045        getConvergenceTest, buildResidual,
1046        petsc.KSPSetConvergenceTest, petsc.KSPConvergedDefault
1047
1048        """
1049        cdef PetscKSPNormType normtype = KSP_NORM_NONE
1050        cdef void* cctx = NULL
1051        cdef PetscBool islsqr = PETSC_FALSE
1052        if converged is not None:
1053            CHKERR(KSPSetConvergenceTest(
1054                    self.ksp, KSP_Converged, NULL, NULL))
1055            if args is None: args = ()
1056            if kargs is None: kargs = {}
1057            self.set_attr('__converged__', (converged, args, kargs))
1058        else:
1059            # this is wrong in general, since different KSP may use
1060            # different convergence tests (like KSPLSQR for example)
1061            # Now we handle LSQR explicitly, but a proper mechanism,
1062            # say KSPGetDefaultConverged would be more appropriate
1063            CHKERR(KSPGetNormType(self.ksp, &normtype))
1064            if normtype != KSP_NORM_NONE:
1065                CHKERR(PetscObjectTypeCompare(<PetscObject>self.ksp,
1066                                              KSPLSQR,  &islsqr))
1067                CHKERR(KSPConvergedDefaultCreate(&cctx))
1068                if not islsqr:
1069                    CHKERR(KSPSetConvergenceTest(self.ksp, KSPConvergedDefault,
1070                                                 cctx, KSPConvergedDefaultDestroy))
1071                else:
1072                    CHKERR(KSPSetConvergenceTest(self.ksp, KSPLSQRConvergedDefault,
1073                                                 cctx, KSPConvergedDefaultDestroy))
1074            else:
1075                CHKERR(KSPSetConvergenceTest(self.ksp, KSPConvergedSkip,
1076                                             NULL, NULL))
1077            self.set_attr('__converged__', None)
1078
1079    def addConvergenceTest(
1080        self,
1081        converged: KSPConvergenceTestFunction,
1082        args: tuple[Any, ...] | None = None,
1083        kargs: dict[str, Any] | None = None,
1084        prepend: bool = False) -> None:
1085        """Add the function to be used to determine convergence.
1086
1087        Logically collective.
1088
1089        Parameters
1090        ----------
1091        converged
1092            Callback which computes the convergence.
1093        args
1094            Positional arguments for callback function.
1095        kargs
1096            Keyword arguments for callback function.
1097        prepend
1098            Whether to prepend this call before the default
1099            convergence test or call it after.
1100
1101        Notes
1102        -----
1103        Cannot be mixed with a call to `setConvergenceTest`.
1104        It can only be called once. If called multiple times, it will
1105        generate an error.
1106
1107        See Also
1108        --------
1109        setTolerances, getConvergenceTest, setConvergenceTest,
1110        petsc.KSPSetConvergenceTest, petsc.KSPConvergedDefault
1111
1112        """
1113        cdef object oconverged = self.get_attr("__converged__")
1114        cdef PetscBool pre = asBool(prepend)
1115        if converged is None: return
1116        if oconverged is not None: raise NotImplementedError("converged callback already set or added")
1117        CHKERR(KSPAddConvergenceTest(self.ksp, KSP_Converged, pre))
1118        if args is None: args = ()
1119        if kargs is None: kargs = {}
1120        self.set_attr('__converged__', (converged, args, kargs))
1121
1122    def getConvergenceTest(self) -> KSPConvergenceTestFunction:
1123        """Return the function to be used to determine convergence.
1124
1125        Logically collective.
1126
1127        See Also
1128        --------
1129        setTolerances, setConvergenceTest, petsc.KSPGetConvergenceTest
1130        petsc.KSPConvergedDefault
1131
1132        """
1133        return self.get_attr('__converged__')
1134
1135    def callConvergenceTest(self, its: int, rnorm: float) -> None:
1136        """Call the convergence test callback.
1137
1138        Collective.
1139
1140        Parameters
1141        ----------
1142        its
1143            Number of iterations.
1144        rnorm
1145            The residual norm.
1146
1147        Notes
1148        -----
1149        This functionality is implemented in petsc4py.
1150
1151        """
1152        cdef PetscInt  ival = asInt(its)
1153        cdef PetscReal rval = asReal(rnorm)
1154        cdef PetscKSPConvergedReason reason = KSP_CONVERGED_ITERATING
1155        CHKERR(KSPConvergenceTestCall(self.ksp, ival, rval, &reason))
1156        return reason
1157
1158    def setConvergenceHistory(
1159        self,
1160        length: int | None = None,
1161        reset: bool = False) -> None:
1162        """Set the array used to hold the residual history.
1163
1164        Not collective.
1165
1166        If set, this array will contain the residual norms computed at
1167        each iteration of the solver.
1168
1169        Parameters
1170        ----------
1171        length
1172            Length of array to store history in.
1173        reset
1174            `True` indicates the history counter is reset to zero for
1175            each new linear solve.
1176
1177        Notes
1178        -----
1179        If ``length`` is not provided or `None` then a default array
1180        of length 10000 is allocated.
1181
1182        If the array is not long enough then once the iterations is
1183        longer than the array length `solve` stops recording the
1184        history.
1185
1186        See Also
1187        --------
1188        getConvergenceHistory, petsc.KSPSetResidualHistory
1189
1190        """
1191        cdef PetscReal *data = NULL
1192        cdef PetscInt   size = 10000
1193        cdef PetscBool flag = PETSC_FALSE
1194        if   length is True:     pass
1195        elif length is not None: size = asInt(length)
1196        if size < 0: size = 10000
1197        if reset: flag = PETSC_TRUE
1198        cdef object hist = oarray_r(empty_r(size), NULL, &data)
1199        self.set_attr('__history__', hist)
1200        CHKERR(KSPSetResidualHistory(self.ksp, data, size, flag))
1201
1202    def getConvergenceHistory(self) -> ArrayReal:
1203        """Return array containing the residual history.
1204
1205        Not collective.
1206
1207        See Also
1208        --------
1209        setConvergenceHistory, petsc.KSPGetResidualHistory
1210
1211        """
1212        cdef const PetscReal *data = NULL
1213        cdef PetscInt   size = 0
1214        CHKERR(KSPGetResidualHistory(self.ksp, &data, &size))
1215        return array_r(size, data)
1216
1217    def logConvergenceHistory(self, rnorm: float) -> None:
1218        """Add residual to convergence history.
1219
1220        Logically collective.
1221
1222        Parameters
1223        ----------
1224        rnorm
1225            Residual norm to be added to convergence history.
1226
1227        """
1228        cdef PetscReal rval = asReal(rnorm)
1229        CHKERR(KSPLogResidualHistory(self.ksp, rval))
1230
1231    # --- monitoring ---
1232
1233    def setMonitor(self,
1234                   monitor: KSPMonitorFunction,
1235                   args: tuple[Any, ...] | None = None,
1236                   kargs: dict[str, Any] | None = None) -> None:
1237        """Set additional function to monitor the residual.
1238
1239        Logically collective.
1240
1241        Set an ADDITIONAL function to be called at every iteration to
1242        monitor the residual/error etc.
1243
1244        Parameters
1245        ----------
1246        monitor
1247            Callback which monitors the convergence.
1248        args
1249            Positional arguments for callback function.
1250        kargs
1251            Keyword arguments for callback function.
1252
1253        Notes
1254        -----
1255        The default is to do nothing. To print the residual, or
1256        preconditioned residual if
1257        ``setNormType(NORM_PRECONDITIONED)`` was called, use
1258        `monitor` as the monitoring routine, with a
1259        `PETSc.Viewer.ASCII` as the context.
1260
1261        Several different monitoring routines may be set by calling
1262        `setMonitor` multiple times; all will be called in the order
1263        in which they were set.
1264
1265        See Also
1266        --------
1267        petsc_options, getMonitor, monitor, monitorCancel, petsc.KSPMonitorSet
1268
1269        """
1270        if monitor is None: return
1271        cdef object monitorlist = self.get_attr('__monitor__')
1272        if monitorlist is None:
1273            monitorlist = []
1274            self.set_attr('__monitor__', monitorlist)
1275            CHKERR(KSPMonitorSet(self.ksp, KSP_Monitor, NULL, NULL))
1276        if args is None: args = ()
1277        if kargs is None: kargs = {}
1278        monitorlist.append((monitor, args, kargs))
1279
1280    def getMonitor(self) -> KSPMonitorFunction:
1281        """Return function used to monitor the residual.
1282
1283        Not collective.
1284
1285        See Also
1286        --------
1287        petsc_options, setMonitor, monitor, monitorCancel
1288        petsc.KSPGetMonitorContext
1289
1290        """
1291        return self.get_attr('__monitor__')
1292
1293    def monitorCancel(self) -> None:
1294        """Clear all monitors for a `KSP` object.
1295
1296        Logically collective.
1297
1298        See Also
1299        --------
1300        petsc_options, getMonitor, setMonitor, monitor, petsc.KSPMonitorCancel
1301
1302        """
1303        CHKERR(KSPMonitorCancel(self.ksp))
1304        self.set_attr('__monitor__', None)
1305
1306    cancelMonitor = monitorCancel
1307
1308    def monitor(self, its: int, rnorm: float) -> None:
1309        """Run the user provided monitor routines, if they exist.
1310
1311        Collective.
1312
1313        Notes
1314        -----
1315        This routine is called by the `KSP` implementations. It does not
1316        typically need to be called by the user.
1317
1318        See Also
1319        --------
1320        setMonitor, petsc.KSPMonitor
1321
1322        """
1323        cdef PetscInt  ival = asInt(its)
1324        cdef PetscReal rval = asReal(rnorm)
1325        CHKERR(KSPMonitor(self.ksp, ival, rval))
1326
1327    # --- customization ---
1328
1329    def setPCSide(self, side: PC.Side) -> None:
1330        """Set the preconditioning side.
1331
1332        Logically collective.
1333
1334        Parameters
1335        ----------
1336        side
1337            The preconditioning side (see `PC.Side`).
1338
1339        Notes
1340        -----
1341        Left preconditioning is used by default for most Krylov methods
1342        except `Type.FGMRES` which only supports right preconditioning.
1343
1344        For methods changing the side of the preconditioner changes the
1345        norm type that is used, see `setNormType`.
1346
1347        Symmetric preconditioning is currently available only for the
1348        `Type.QCG` method. Note, however, that symmetric preconditioning
1349        can be emulated by using either right or left preconditioning
1350        and a pre or post processing step.
1351
1352        Setting the PC side often affects the default norm type. See
1353        `setNormType` for details.
1354
1355        See Also
1356        --------
1357        PC.Side, petsc_options, getPCSide, setNormType, getNormType
1358        petsc.KSPSetPCSide
1359
1360        """
1361        CHKERR(KSPSetPCSide(self.ksp, side))
1362
1363    def getPCSide(self) -> PC.Side:
1364        """Return the preconditioning side.
1365
1366        Not collective.
1367
1368        See Also
1369        --------
1370        petsc_options, setPCSide, setNormType, getNormType, petsc.KSPGetPCSide
1371
1372        """
1373        cdef PetscPCSide side = PC_LEFT
1374        CHKERR(KSPGetPCSide(self.ksp, &side))
1375        return side
1376
1377    def setNormType(self, normtype: NormType) -> None:
1378        """Set the norm that is used for convergence testing.
1379
1380        Logically collective.
1381
1382        Parameters
1383        ----------
1384        normtype
1385            The norm type to use (see `NormType`).
1386
1387        Notes
1388        -----
1389        Not all combinations of preconditioner side (see
1390        `setPCSide`) and norm type are supported by all Krylov
1391        methods. If only one is set, PETSc tries to automatically
1392        change the other to find a compatible pair. If no such
1393        combination is supported, PETSc will generate an error.
1394
1395        See Also
1396        --------
1397        NormType, petsc_options, setUp, solve, destroy, setPCSide, getPCSide
1398        NormType, petsc.KSPSetNormType, petsc.KSPConvergedSkip
1399        petsc.KSPSetCheckNormIteration
1400
1401        """
1402        CHKERR(KSPSetNormType(self.ksp, normtype))
1403
1404    def getNormType(self) -> NormType:
1405        """Return the norm that is used for convergence testing.
1406
1407        Not collective.
1408
1409        See Also
1410        --------
1411        NormType, setNormType, petsc.KSPGetNormType, petsc.KSPConvergedSkip
1412
1413        """
1414        cdef PetscKSPNormType normtype = KSP_NORM_NONE
1415        CHKERR(KSPGetNormType(self.ksp, &normtype))
1416        return normtype
1417
1418    def setComputeEigenvalues(self, flag: bool) -> None:
1419        """Set a flag to compute eigenvalues.
1420
1421        Logically collective.
1422
1423        Set a flag so that the extreme eigenvalues values will be
1424        calculated via a Lanczos or Arnoldi process as the linear
1425        system is solved.
1426
1427        Parameters
1428        ----------
1429        flag
1430            Boolean whether to compute eigenvalues (or not).
1431
1432        Notes
1433        -----
1434        Currently this option is not valid for all iterative methods.
1435
1436        See Also
1437        --------
1438        getComputeEigenvalues, petsc.KSPSetComputeEigenvalues
1439
1440        """
1441        cdef PetscBool compute = asBool(flag)
1442        CHKERR(KSPSetComputeEigenvalues(self.ksp, compute))
1443
1444    def getComputeEigenvalues(self) -> bool:
1445        """Return flag indicating whether eigenvalues will be calculated.
1446
1447        Not collective.
1448
1449        Return the flag indicating that the extreme eigenvalues values
1450        will be calculated via a Lanczos or Arnoldi process as the
1451        linear system is solved.
1452
1453        See Also
1454        --------
1455        setComputeEigenvalues, petsc.KSPSetComputeEigenvalues
1456
1457        """
1458        cdef PetscBool flag = PETSC_FALSE
1459        CHKERR(KSPGetComputeEigenvalues(self.ksp, &flag))
1460        return toBool(flag)
1461
1462    def setComputeSingularValues(self, flag: bool) -> None:
1463        """Set flag to calculate singular values.
1464
1465        Logically collective.
1466
1467        Set a flag so that the extreme singular values will be
1468        calculated via a Lanczos or Arnoldi process as the linear
1469        system is solved.
1470
1471        Parameters
1472        ----------
1473        flag
1474            Boolean whether to compute singular values (or not).
1475
1476        Notes
1477        -----
1478        Currently this option is not valid for all iterative methods.
1479
1480        See Also
1481        --------
1482        getComputeSingularValues, petsc.KSPSetComputeSingularValues
1483
1484        """
1485        cdef PetscBool compute = asBool(flag)
1486        CHKERR(KSPSetComputeSingularValues(self.ksp, compute))
1487
1488    def getComputeSingularValues(self) -> bool:
1489        """Return flag indicating whether singular values will be calculated.
1490
1491        Not collective.
1492
1493        Return the flag indicating whether the extreme singular values
1494        will be calculated via a Lanczos or Arnoldi process as the
1495        linear system is solved.
1496
1497        See Also
1498        --------
1499        setComputeSingularValues, petsc.KSPGetComputeSingularValues
1500
1501        """
1502        cdef PetscBool flag = PETSC_FALSE
1503        CHKERR(KSPGetComputeSingularValues(self.ksp, &flag))
1504        return toBool(flag)
1505
1506    # --- initial guess ---
1507
1508    def setInitialGuessNonzero(self, flag: bool) -> None:
1509        """Tell the iterative solver that the initial guess is nonzero.
1510
1511        Logically collective.
1512
1513        Otherwise KSP assumes the initial guess is to be zero (and thus
1514        zeros it out before solving).
1515
1516        Parameters
1517        ----------
1518        flag
1519            `True` indicates the guess is non-zero, `False`
1520            indicates the guess is zero.
1521
1522        See Also
1523        --------
1524        petsc.KSPSetInitialGuessNonzero
1525
1526        """
1527        cdef PetscBool guess_nonzero = asBool(flag)
1528        CHKERR(KSPSetInitialGuessNonzero(self.ksp, guess_nonzero))
1529
1530    def getInitialGuessNonzero(self) -> bool:
1531        """Determine whether the KSP solver uses a zero initial guess.
1532
1533        Not collective.
1534
1535        See Also
1536        --------
1537        petsc.KSPGetInitialGuessNonzero
1538
1539        """
1540        cdef PetscBool flag = PETSC_FALSE
1541        CHKERR(KSPGetInitialGuessNonzero(self.ksp, &flag))
1542        return toBool(flag)
1543
1544    def setInitialGuessKnoll(self, flag: bool) -> None:
1545        """Tell solver to use `PC.apply` to compute the initial guess.
1546
1547        Logically collective.
1548
1549        This is the Knoll trick.
1550
1551        Parameters
1552        ----------
1553        flag
1554            `True` uses Knoll trick.
1555
1556        See Also
1557        --------
1558        petsc.KSPSetInitialGuessKnoll
1559
1560        """
1561        cdef PetscBool guess_knoll = asBool(flag)
1562        CHKERR(KSPSetInitialGuessKnoll(self.ksp, guess_knoll))
1563
1564    def getInitialGuessKnoll(self) -> bool:
1565        """Determine whether the KSP solver is using the Knoll trick.
1566
1567        Not collective.
1568
1569        This uses the Knoll trick; using `PC.apply` to compute the
1570        initial guess.
1571
1572        See Also
1573        --------
1574        petsc.KSPGetInitialGuessKnoll
1575
1576        """
1577        cdef PetscBool flag = PETSC_FALSE
1578        CHKERR(KSPGetInitialGuessKnoll(self.ksp, &flag))
1579        return toBool(flag)
1580
1581    def setUseFischerGuess(self, model: int, size: int) -> None:
1582        """Use the Paul Fischer algorithm to compute initial guesses.
1583
1584        Logically collective.
1585
1586        Use the Paul Fischer algorithm or its variants to compute
1587        initial guesses for a set of solves with related right hand
1588        sides.
1589
1590        Parameters
1591        ----------
1592        model
1593            Use model ``1``, model ``2``, model ``3``, any other number
1594            to turn it off.
1595        size
1596            Size of subspace used to generate initial guess.
1597
1598        See Also
1599        --------
1600        petsc.KSPSetUseFischerGuess
1601
1602        """
1603        cdef PetscInt ival1 = asInt(model)
1604        cdef PetscInt ival2 = asInt(size)
1605        CHKERR(KSPSetUseFischerGuess(self.ksp, ival1, ival2))
1606
1607    # --- solving ---
1608
1609    def setUp(self) -> None:
1610        """Set up internal data structures for an iterative solver.
1611
1612        Collective.
1613
1614        See Also
1615        --------
1616        petsc.KSPSetUp
1617
1618        """
1619        CHKERR(KSPSetUp(self.ksp))
1620
1621    def reset(self) -> None:
1622        """Resets a KSP context.
1623
1624        Collective.
1625
1626        Resets a KSP context to the ``kspsetupcalled = 0`` state and
1627        removes any allocated Vecs and Mats.
1628
1629        See Also
1630        --------
1631        petsc.KSPReset
1632
1633        """
1634        CHKERR(KSPReset(self.ksp))
1635
1636    def setUpOnBlocks(self) -> None:
1637        """Set up the preconditioner for each block in a block method.
1638
1639        Collective.
1640
1641        Methods include: block Jacobi, block Gauss-Seidel, and
1642        overlapping Schwarz methods.
1643
1644        See Also
1645        --------
1646        petsc.KSPSetUpOnBlocks
1647
1648        """
1649        CHKERR(KSPSetUpOnBlocks(self.ksp))
1650
1651    def setPreSolve(
1652        self,
1653        presolve: KSPPreSolveFunction | None,
1654        args: tuple[Any, ...] | None = None,
1655        kargs: dict[str, Any] | None = None) -> None:
1656        """Set the function that is called at the beginning of each `KSP.solve`.
1657
1658        Logically collective.
1659
1660        Parameters
1661        ----------
1662        presolve
1663            The callback function.
1664        args
1665            Positional arguments for the callback function.
1666        kargs
1667            Keyword arguments for the callback function.
1668
1669        See Also
1670        --------
1671        solve, petsc.KSPSetPreSolve, petsc.KSPSetPostSolve
1672
1673        """
1674        if presolve is not None:
1675            if args is None: args = ()
1676            if kargs is None: kargs = {}
1677            context = (presolve, args, kargs)
1678            self.set_attr('__presolve__', context)
1679            CHKERR(KSPSetPreSolve(self.ksp, KSP_PreSolve, <void*>context))
1680        else:
1681            self.set_attr('__presolve__', None)
1682            CHKERR(KSPSetPreSolve(self.ksp, NULL, NULL))
1683
1684    def setPostSolve(
1685        self,
1686        postsolve: KSPPostSolveFunction | None,
1687        args: tuple[Any, ...] | None = None,
1688        kargs: dict[str, Any] | None = None) -> None:
1689        """Set the function that is called at the end of each `KSP.solve`.
1690
1691        Logically collective.
1692
1693        Parameters
1694        ----------
1695        postsolve
1696            The callback function.
1697        args
1698            Positional arguments for the callback function.
1699        kargs
1700            Keyword arguments for the callback function.
1701
1702        See Also
1703        --------
1704        solve, petsc.KSPSetPreSolve, petsc.KSPSetPostSolve
1705
1706        """
1707        if postsolve is not None:
1708            if args is None: args = ()
1709            if kargs is None: kargs = {}
1710            context = (postsolve, args, kargs)
1711            self.set_attr('__postsolve__', context)
1712            CHKERR(KSPSetPostSolve(self.ksp, KSP_PostSolve, <void*>context))
1713        else:
1714            self.set_attr('__postsolve__', None)
1715            CHKERR(KSPSetPostSolve(self.ksp, NULL, NULL))
1716
1717    def solve(self, Vec b, Vec x) -> None:
1718        """Solve the linear system.
1719
1720        Collective.
1721
1722        Parameters
1723        ----------
1724        b
1725            Right hand side vector.
1726        x
1727            Solution vector.
1728
1729        Notes
1730        -----
1731        If one uses `setDM` then ``x`` or ``b`` need not be passed. Use
1732        `getSolution` to access the solution in this case.
1733
1734        The operator is specified with `setOperators`.
1735
1736        `solve` will normally return without generating an error
1737        regardless of whether the linear system was solved or if
1738        constructing the preconditioner failed. Call
1739        `getConvergedReason` to determine if the solver converged or
1740        failed and why. The option ``-ksp_error_if_not_converged`` or
1741        function `setErrorIfNotConverged` will cause `solve` to error
1742        as soon as an error occurs in the linear solver. In inner
1743        solves, ``DIVERGED_MAX_IT`` is not treated as an error
1744        because when using nested solvers it may be fine that inner
1745        solvers in the preconditioner do not converge during the
1746        solution process.
1747
1748        The number of iterations can be obtained from `its`.
1749
1750        If you provide a matrix that has a `Mat.setNullSpace` and
1751        `Mat.setTransposeNullSpace` this will use that information to
1752        solve singular systems in the least squares sense with a norm
1753        minimizing solution.
1754
1755        Ax = b where b = bₚ + bₜ where bₜ is not in the range of A
1756        (and hence by the fundamental theorem of linear algebra is in
1757        the nullspace(Aᵀ), see `Mat.setNullSpace`.
1758
1759        KSP first removes bₜ producing the linear system Ax = bₚ (which
1760        has multiple solutions) and solves this to find the ∥x∥
1761        minimizing solution (and hence it finds the solution x
1762        orthogonal to the nullspace(A). The algorithm is simply in each
1763        iteration of the Krylov method we remove the nullspace(A) from
1764        the search direction thus the solution which is a linear
1765        combination of the search directions has no component in the
1766        nullspace(A).
1767
1768        We recommend always using `Type.GMRES` for such singular
1769        systems. If nullspace(A) = nullspace(Aᵀ) (note symmetric
1770        matrices always satisfy this property) then both left and right
1771        preconditioning will work If nullspace(A) != nullspace(Aᵀ) then
1772        left preconditioning will work but right preconditioning may
1773        not work (or it may).
1774
1775        If using a direct method (e.g., via the KSP solver
1776        `Type.PREONLY` and a preconditioner such as `PC.Type.LU` or
1777        `PC.Type.ILU`, then its=1. See `setTolerances` for more details.
1778
1779        **Understanding Convergence**
1780
1781        The routines `setMonitor` and `computeEigenvalues` provide
1782        information on additional options to monitor convergence and
1783        print eigenvalue information.
1784
1785        See Also
1786        --------
1787        create, setUp, destroy, setTolerances, is_converged, solveTranspose, its
1788        Mat.setNullSpace, Mat.setTransposeNullSpace, Type,
1789        setErrorIfNotConverged petsc_options, petsc.KSPSolve
1790
1791        """
1792        cdef PetscVec b_vec = NULL
1793        cdef PetscVec x_vec = NULL
1794        if b is not None: b_vec = b.vec
1795        if x is not None: x_vec = x.vec
1796        CHKERR(KSPSolve(self.ksp, b_vec, x_vec))
1797
1798    def solveTranspose(self, Vec b, Vec x) -> None:
1799        """Solve the transpose of a linear system.
1800
1801        Collective.
1802
1803        Parameters
1804        ----------
1805        b
1806            Right hand side vector.
1807        x
1808            Solution vector.
1809
1810        Notes
1811        -----
1812        For complex numbers this solve the non-Hermitian transpose
1813        system.
1814
1815        See Also
1816        --------
1817        solve, petsc.KSPSolveTranspose
1818
1819        """
1820        CHKERR(KSPSolveTranspose(self.ksp, b.vec, x.vec))
1821
1822    def matSolve(self, Mat B, Mat X) -> None:
1823        """Solve a linear system with multiple right-hand sides.
1824
1825        Collective.
1826
1827        These are stored as a `Mat.Type.DENSE`. Unlike `solve`,
1828        ``B`` and ``X`` must be different matrices.
1829
1830        Parameters
1831        ----------
1832        B
1833            Block of right-hand sides.
1834        X
1835            Block of solutions.
1836
1837        See Also
1838        --------
1839        solve, petsc.KSPMatSolve
1840
1841        """
1842        CHKERR(KSPMatSolve(self.ksp, B.mat, X.mat))
1843
1844    def matSolveTranspose(self, Mat B, Mat X) -> None:
1845        """Solve the transpose of a linear system with multiple RHS.
1846
1847        Collective.
1848
1849        Parameters
1850        ----------
1851        B
1852            Block of right-hand sides.
1853        X
1854            Block of solutions.
1855
1856        See Also
1857        --------
1858        solveTranspose, petsc.KSPMatSolve
1859
1860        """
1861        CHKERR(KSPMatSolveTranspose(self.ksp, B.mat, X.mat))
1862
1863    def setIterationNumber(self, its: int) -> None:
1864        """Use `its` property."""
1865        cdef PetscInt ival = asInt(its)
1866        CHKERR(KSPSetIterationNumber(self.ksp, ival))
1867
1868    def getIterationNumber(self) -> int:
1869        """Use `its` property."""
1870        cdef PetscInt ival = 0
1871        CHKERR(KSPGetIterationNumber(self.ksp, &ival))
1872        return toInt(ival)
1873
1874    def setResidualNorm(self, rnorm: float) -> None:
1875        """Use `norm` property."""
1876        cdef PetscReal rval = asReal(rnorm)
1877        CHKERR(KSPSetResidualNorm(self.ksp, rval))
1878
1879    def getResidualNorm(self) -> float:
1880        """Use `norm` property."""
1881        cdef PetscReal rval = 0
1882        CHKERR(KSPGetResidualNorm(self.ksp, &rval))
1883        return toReal(rval)
1884
1885    def setConvergedReason(self, reason: KSP.ConvergedReason) -> None:
1886        """Use `reason` property."""
1887        cdef PetscKSPConvergedReason val = reason
1888        CHKERR(KSPSetConvergedReason(self.ksp, val))
1889
1890    def getConvergedReason(self) -> KSP.ConvergedReason:
1891        """Use `reason` property."""
1892        cdef PetscKSPConvergedReason reason = KSP_CONVERGED_ITERATING
1893        CHKERR(KSPGetConvergedReason(self.ksp, &reason))
1894        return reason
1895
1896    def getCGObjectiveValue(self) -> float:
1897        """Return the CG objective function value.
1898
1899        Not collective.
1900
1901        See Also
1902        --------
1903        petsc.KSPCGGetObjFcn
1904
1905        """
1906        cdef PetscReal cval = 0
1907        CHKERR(KSPCGGetObjFcn(self.ksp, &cval))
1908        return cval
1909
1910    def setHPDDMType(self, hpddm_type: HPDDMType) -> None:
1911        """Set the Krylov solver type.
1912
1913        Collective.
1914
1915        Parameters
1916        ----------
1917        hpddm_type
1918            The type of Krylov solver to use.
1919
1920        See Also
1921        --------
1922        petsc.KSPHPDDMSetType
1923
1924        """
1925        cdef PetscKSPHPDDMType ctype = hpddm_type
1926        CHKERR(KSPHPDDMSetType(self.ksp, ctype))
1927
1928    def getHPDDMType(self) -> HPDDMType:
1929        """Return the Krylov solver type.
1930
1931        Not collective.
1932
1933        See Also
1934        --------
1935        petsc.KSPHPDDMGetType
1936
1937        """
1938        cdef PetscKSPHPDDMType cval = KSP_HPDDM_TYPE_GMRES
1939        CHKERR(KSPHPDDMGetType(self.ksp, &cval))
1940        return cval
1941
1942    def setErrorIfNotConverged(self, flag: bool) -> None:
1943        """Cause `solve` to generate an error if not converged.
1944
1945        Logically collective.
1946
1947        Parameters
1948        ----------
1949        flag
1950            `True` enables this behavior.
1951
1952        See Also
1953        --------
1954        petsc.KSPSetErrorIfNotConverged
1955
1956        """
1957        cdef PetscBool ernc = asBool(flag)
1958        CHKERR(KSPSetErrorIfNotConverged(self.ksp, ernc))
1959
1960    def getErrorIfNotConverged(self) -> bool:
1961        """Return the flag indicating the solver will error if divergent.
1962
1963        Not collective.
1964
1965        See Also
1966        --------
1967        petsc.KSPGetErrorIfNotConverged
1968
1969        """
1970        cdef PetscBool flag = PETSC_FALSE
1971        CHKERR(KSPGetErrorIfNotConverged(self.ksp, &flag))
1972        return toBool(flag)
1973
1974    def getRhs(self) -> Vec:
1975        """Return the right-hand side vector for the linear system.
1976
1977        Not collective.
1978
1979        See Also
1980        --------
1981        petsc.KSPGetRhs
1982
1983        """
1984        cdef Vec vec = Vec()
1985        CHKERR(KSPGetRhs(self.ksp, &vec.vec))
1986        CHKERR(PetscINCREF(vec.obj))
1987        return vec
1988
1989    def getSolution(self) -> Vec:
1990        """Return the solution for the linear system to be solved.
1991
1992        Not collective.
1993
1994        Note that this may not be the solution that is stored during
1995        the iterative process.
1996
1997        See Also
1998        --------
1999        petsc.KSPGetSolution
2000
2001        """
2002        cdef Vec vec = Vec()
2003        CHKERR(KSPGetSolution(self.ksp, &vec.vec))
2004        CHKERR(PetscINCREF(vec.obj))
2005        return vec
2006
2007    def getWorkVecs(
2008        self,
2009        right: int | None = None,
2010        left: int | None = None) -> tuple[list[Vec], list[Vec]] | list[Vec] | None:
2011        """Create working vectors.
2012
2013        Collective.
2014
2015        Parameters
2016        ----------
2017        right
2018            Number of right hand vectors to allocate.
2019        left
2020            Number of left hand vectors to allocate.
2021
2022        Returns
2023        -------
2024        R : list of Vec
2025            List of correctly allocated right hand vectors.
2026        L : list of Vec
2027            List of correctly allocated left hand vectors.
2028
2029        """
2030        cdef bint R = right is not None
2031        cdef bint L = left  is not None
2032        cdef PetscInt i=0, nr=0, nl=0
2033        cdef PetscVec *vr=NULL, *vl=NULL
2034        if R: nr = asInt(right)
2035        if L: nl = asInt(left)
2036        cdef object vecsr = [] if R else None
2037        cdef object vecsl = [] if L else None
2038        CHKERR(KSPCreateVecs(self.ksp, nr, &vr, nl, &vr))
2039        try:
2040            for i from 0 <= i < nr:
2041                vecsr.append(ref_Vec(vr[i]))
2042            for i from 0 <= i < nl:
2043                vecsl.append(ref_Vec(vl[i]))
2044        finally:
2045            if nr > 0 and vr != NULL:
2046                VecDestroyVecs(nr, &vr) # XXX errors?
2047            if nl > 0 and vl !=NULL:
2048                VecDestroyVecs(nl, &vl) # XXX errors?
2049        #
2050        if R and L: return (vecsr, vecsl)
2051        elif R:     return vecsr
2052        elif L:     return vecsl
2053        else:       return None
2054
2055    def buildSolution(self, Vec x=None) -> Vec:
2056        """Return the solution vector.
2057
2058        Collective.
2059
2060        Parameters
2061        ----------
2062        x
2063            Optional vector to store the solution.
2064
2065        See Also
2066        --------
2067        buildResidual, petsc.KSPBuildSolution
2068
2069        """
2070        if x is None: x = Vec()
2071        if x.vec == NULL:
2072            CHKERR(KSPGetSolution(self.ksp, &x.vec))
2073            CHKERR(VecDuplicate(x.vec, &x.vec))
2074        CHKERR(KSPBuildSolution(self.ksp, x.vec, NULL))
2075        return x
2076
2077    def buildResidual(self, Vec r=None) -> Vec:
2078        """Return the residual of the linear system.
2079
2080        Collective.
2081
2082        Parameters
2083        ----------
2084        r
2085            Optional vector to use for the result.
2086
2087        See Also
2088        --------
2089        buildSolution, petsc.KSPBuildResidual
2090
2091        """
2092        if r is None: r = Vec()
2093        if r.vec == NULL:
2094            CHKERR(KSPGetRhs(self.ksp, &r.vec))
2095            CHKERR(VecDuplicate(r.vec, &r.vec))
2096        CHKERR(KSPBuildResidual(self.ksp , NULL, r.vec, &r.vec))
2097        return r
2098
2099    def computeEigenvalues(self) -> ArrayComplex:
2100        """Compute the extreme eigenvalues for the preconditioned operator.
2101
2102        Not collective.
2103
2104        See Also
2105        --------
2106        petsc.KSPComputeEigenvalues
2107
2108        """
2109        cdef PetscInt its = 0
2110        cdef PetscInt neig = 0
2111        cdef PetscReal *rdata = NULL
2112        cdef PetscReal *idata = NULL
2113        CHKERR(KSPGetIterationNumber(self.ksp, &its))
2114        cdef ndarray r = oarray_r(empty_r(its), NULL, &rdata)
2115        cdef ndarray i = oarray_r(empty_r(its), NULL, &idata)
2116        CHKERR(KSPComputeEigenvalues(self.ksp, its, rdata, idata, &neig))
2117        eigen = empty_c(neig)
2118        eigen.real = r[:neig]
2119        eigen.imag = i[:neig]
2120        return eigen
2121
2122    def computeExtremeSingularValues(self) -> tuple[float, float]:
2123        """Compute the extreme singular values for the preconditioned operator.
2124
2125        Collective.
2126
2127        Returns
2128        -------
2129        smax : float
2130            The maximum singular value.
2131        smin : float
2132            The minimum singular value.
2133
2134        See Also
2135        --------
2136        petsc.KSPComputeExtremeSingularValues
2137
2138        """
2139        cdef PetscReal smax = 0
2140        cdef PetscReal smin = 0
2141        CHKERR(KSPComputeExtremeSingularValues(self.ksp, &smax, &smin))
2142        return toReal(smax), toReal(smin)
2143
2144    # --- GMRES ---
2145
2146    def setGMRESRestart(self, restart: int) -> None:
2147        """Set number of iterations at which KSP restarts.
2148
2149        Logically collective.
2150
2151        Suitable KSPs are: KSPGMRES, KSPFGMRES and KSPLGMRES.
2152
2153        Parameters
2154        ----------
2155        restart
2156            Integer restart value.
2157
2158        See Also
2159        --------
2160        petsc.KSPGMRESSetRestart
2161
2162        """
2163        cdef PetscInt ival = asInt(restart)
2164        CHKERR(KSPGMRESSetRestart(self.ksp, ival))
2165
2166    # --- Python ---
2167
2168    def createPython(
2169        self,
2170        context: Any = None,
2171        comm: Comm | None = None) -> Self:
2172        """Create a linear solver of Python type.
2173
2174        Collective.
2175
2176        Parameters
2177        ----------
2178        context
2179            An instance of the Python class implementing the required
2180            methods.
2181        comm
2182            MPI communicator, defaults to `Sys.getDefaultComm`.
2183
2184        See Also
2185        --------
2186        petsc_python_ksp, setType, setPythonContext, Type.PYTHON
2187
2188        """
2189        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_DEFAULT)
2190        cdef PetscKSP newksp = NULL
2191        CHKERR(KSPCreate(ccomm, &newksp))
2192        CHKERR(PetscCLEAR(self.obj)); self.ksp = newksp
2193        CHKERR(KSPSetType(self.ksp, KSPPYTHON))
2194        CHKERR(KSPPythonSetContext(self.ksp, <void*>context))
2195        return self
2196
2197    def setPythonContext(self, context: Any | None = None) -> None:
2198        """Set the instance of the class implementing Python methods.
2199
2200        Not collective.
2201
2202        See Also
2203        --------
2204        petsc_python_ksp, getPythonContext
2205
2206        """
2207        CHKERR(KSPPythonSetContext(self.ksp, <void*>context))
2208
2209    def getPythonContext(self) -> Any:
2210        """Return the instance of the class implementing Python methods.
2211
2212        Not collective.
2213
2214        See Also
2215        --------
2216        petsc_python_ksp, setPythonContext
2217
2218        """
2219        cdef void *context = NULL
2220        CHKERR(KSPPythonGetContext(self.ksp, &context))
2221        if context == NULL: return None
2222        else: return <object> context
2223
2224    def setPythonType(self, py_type: str) -> None:
2225        """Set the fully qualified Python name of the class to be used.
2226
2227        Collective.
2228
2229        See Also
2230        --------
2231        petsc_python_ksp, setPythonContext, getPythonType
2232        petsc.KSPPythonSetType
2233
2234        """
2235        cdef const char *cval = NULL
2236        py_type = str2bytes(py_type, &cval)
2237        CHKERR(KSPPythonSetType(self.ksp, cval))
2238
2239    def getPythonType(self) -> str:
2240        """Return the fully qualified Python name of the class used by the solver.
2241
2242        Not collective.
2243
2244        See Also
2245        --------
2246        petsc_python_ksp, setPythonContext, setPythonType
2247        petsc.KSPPythonGetType
2248
2249        """
2250        cdef const char *cval = NULL
2251        CHKERR(KSPPythonGetType(self.ksp, &cval))
2252        return bytes2str(cval)
2253
2254    # --- application context ---
2255
2256    property appctx:
2257        """The solver application context."""
2258        def __get__(self) -> Any:
2259            return self.getAppCtx()
2260
2261        def __set__(self, value):
2262            self.setAppCtx(value)
2263
2264    # --- discretization space ---
2265
2266    property dm:
2267        """The solver `DM`."""
2268        def __get__(self) -> DM:
2269            return self.getDM()
2270
2271        def __set__(self, value):
2272            self.setDM(value)
2273
2274    # --- vectors ---
2275
2276    property vec_sol:
2277        """The solution vector."""
2278        def __get__(self) -> Vec:
2279            return self.getSolution()
2280
2281    property vec_rhs:
2282        """The right-hand side vector."""
2283        def __get__(self) -> Vec:
2284            return self.getRhs()
2285
2286    # --- operators ---
2287
2288    property mat_op:
2289        """The system matrix operator."""
2290        def __get__(self) -> Mat:
2291            return self.getOperators()[0]
2292
2293    property mat_pc:
2294        """The preconditioner operator."""
2295        def __get__(self) -> Mat:
2296            return self.getOperators()[1]
2297
2298    # --- initial guess ---
2299
2300    property guess_nonzero:
2301        """Whether guess is non-zero."""
2302        def __get__(self) -> bool:
2303            return self.getInitialGuessNonzero()
2304
2305        def __set__(self, value):
2306            self.setInitialGuessNonzero(value)
2307
2308    property guess_knoll:
2309        """Whether solver uses Knoll trick."""
2310        def __get__(self) -> bool:
2311            return self.getInitialGuessKnoll()
2312
2313        def __set__(self, value):
2314            self.setInitialGuessKnoll(value)
2315
2316    # --- preconditioner ---
2317
2318    property pc:
2319        """The `PC` of the solver."""
2320        def __get__(self) -> PC:
2321            return self.getPC()
2322
2323    property pc_side:
2324        """The side on which preconditioning is performed."""
2325        def __get__(self) -> PC.Side:
2326            return self.getPCSide()
2327
2328        def __set__(self, value):
2329            self.setPCSide(value)
2330
2331    property norm_type:
2332        """The norm used by the solver."""
2333        def __get__(self) -> NormType:
2334            return self.getNormType()
2335
2336        def __set__(self, value):
2337            self.setNormType(value)
2338
2339    # --- tolerances ---
2340
2341    property rtol:
2342        """The relative tolerance of the solver."""
2343        def __get__(self) -> float:
2344            return self.getTolerances()[0]
2345
2346        def __set__(self, value):
2347            self.setTolerances(rtol=value)
2348
2349    property atol:
2350        """The absolute tolerance of the solver."""
2351        def __get__(self) -> float:
2352            return self.getTolerances()[1]
2353
2354        def __set__(self, value):
2355            self.setTolerances(atol=value)
2356
2357    property divtol:
2358        """The divergence tolerance of the solver."""
2359        def __get__(self) -> float:
2360            return self.getTolerances()[2]
2361
2362        def __set__(self, value):
2363            self.setTolerances(divtol=value)
2364
2365    property max_it:
2366        """The maximum number of iteration the solver may take."""
2367        def __get__(self) -> int:
2368            return self.getTolerances()[3]
2369
2370        def __set__(self, value):
2371            self.setTolerances(max_it=value)
2372
2373    # --- iteration ---
2374
2375    property its:
2376        """The current number of iterations the solver has taken."""
2377        def __get__(self) -> int:
2378            return self.getIterationNumber()
2379
2380        def __set__(self, value):
2381            self.setIterationNumber(value)
2382
2383    property norm:
2384        """The norm of the residual at the current iteration."""
2385        def __get__(self) -> float:
2386            return self.getResidualNorm()
2387
2388        def __set__(self, value):
2389            self.setResidualNorm(value)
2390
2391    property history:
2392        """The convergence history of the solver."""
2393        def __get__(self) -> ndarray:
2394            return self.getConvergenceHistory()
2395
2396    # --- convergence ---
2397
2398    property reason:
2399        """The converged reason."""
2400        def __get__(self) -> KSP.ConvergedReason:
2401            return self.getConvergedReason()
2402
2403        def __set__(self, value):
2404            self.setConvergedReason(value)
2405
2406    property is_iterating:
2407        """Boolean indicating if the solver has not converged yet."""
2408        def __get__(self) -> bool:
2409            return self.reason == 0
2410
2411    property is_converged:
2412        """Boolean indicating if the solver has converged."""
2413        def __get__(self) -> bool:
2414            return self.reason > 0
2415
2416    property is_diverged:
2417        """Boolean indicating if the solver has failed."""
2418        def __get__(self) -> bool:
2419            return self.reason < 0
2420
2421# --------------------------------------------------------------------
2422
2423del KSPType
2424del KSPNormType
2425del KSPConvergedReason
2426del KSPHPDDMType
2427
2428# --------------------------------------------------------------------
2429