xref: /petsc/src/binding/petsc4py/src/petsc4py/PETSc/SNES.pyx (revision bd89dbf26d8a5efecb980364933175da61864cd7)
1# --------------------------------------------------------------------
2
3class SNESType(object):
4    """SNES solver type.
5
6    See Also
7    --------
8    petsc.SNESType
9
10    """
11    NEWTONLS         = S_(SNESNEWTONLS)
12    NEWTONTR         = S_(SNESNEWTONTR)
13    NEWTONAL         = S_(SNESNEWTONAL)
14    PYTHON           = S_(SNESPYTHON)
15    NRICHARDSON      = S_(SNESNRICHARDSON)
16    KSPONLY          = S_(SNESKSPONLY)
17    KSPTRANSPOSEONLY = S_(SNESKSPTRANSPOSEONLY)
18    VINEWTONRSLS     = S_(SNESVINEWTONRSLS)
19    VINEWTONSSLS     = S_(SNESVINEWTONSSLS)
20    NGMRES           = S_(SNESNGMRES)
21    QN               = S_(SNESQN)
22    SHELL            = S_(SNESSHELL)
23    NGS              = S_(SNESNGS)
24    NCG              = S_(SNESNCG)
25    FAS              = S_(SNESFAS)
26    MS               = S_(SNESMS)
27    NASM             = S_(SNESNASM)
28    ANDERSON         = S_(SNESANDERSON)
29    ASPIN            = S_(SNESASPIN)
30    COMPOSITE        = S_(SNESCOMPOSITE)
31    PATCH            = S_(SNESPATCH)
32
33
34class SNESNormSchedule(object):
35    """SNES norm schedule.
36
37    See Also
38    --------
39    petsc.SNESNormSchedule
40
41    """
42    # native
43    NORM_DEFAULT            = SNES_NORM_DEFAULT
44    NORM_NONE               = SNES_NORM_NONE
45    NORM_ALWAYS             = SNES_NORM_ALWAYS
46    NORM_INITIAL_ONLY       = SNES_NORM_INITIAL_ONLY
47    NORM_FINAL_ONLY         = SNES_NORM_FINAL_ONLY
48    NORM_INITIAL_FINAL_ONLY = SNES_NORM_INITIAL_FINAL_ONLY
49    # aliases
50    DEFAULT            = NORM_DEFAULT
51    NONE               = NORM_NONE
52    ALWAYS             = NORM_ALWAYS
53    INITIAL_ONLY       = NORM_INITIAL_ONLY
54    FINAL_ONLY         = NORM_FINAL_ONLY
55    INITIAL_FINAL_ONLY = NORM_INITIAL_FINAL_ONLY
56
57
58# FIXME Missing reference petsc.SNESConvergedReason
59class SNESConvergedReason(object):
60    """SNES solver termination reason.
61
62    See Also
63    --------
64    petsc.SNESGetConvergedReason
65
66    """
67    # iterating
68    CONVERGED_ITERATING         = SNES_CONVERGED_ITERATING
69    ITERATING                   = SNES_CONVERGED_ITERATING
70    # converged
71    CONVERGED_FNORM_ABS         = SNES_CONVERGED_FNORM_ABS
72    CONVERGED_FNORM_RELATIVE    = SNES_CONVERGED_FNORM_RELATIVE
73    CONVERGED_SNORM_RELATIVE    = SNES_CONVERGED_SNORM_RELATIVE
74    CONVERGED_ITS               = SNES_CONVERGED_ITS
75    # diverged
76    DIVERGED_FUNCTION_DOMAIN    = SNES_DIVERGED_FUNCTION_DOMAIN
77    DIVERGED_FUNCTION_NANORINF  = SNES_DIVERGED_FUNCTION_NANORINF
78    DIVERGED_OBJECTIVE_DOMAIN   = SNES_DIVERGED_OBJECTIVE_DOMAIN
79    DIVERGED_OBJECTIVE_NANORINF = SNES_DIVERGED_OBJECTIVE_NANORINF
80    DIVERGED_JACOBIAN_DOMAIN    = SNES_DIVERGED_JACOBIAN_DOMAIN
81    DIVERGED_FUNCTION_COUNT     = SNES_DIVERGED_FUNCTION_COUNT
82    DIVERGED_LINEAR_SOLVE       = SNES_DIVERGED_LINEAR_SOLVE
83    DIVERGED_MAX_IT             = SNES_DIVERGED_MAX_IT
84    DIVERGED_LINE_SEARCH        = SNES_DIVERGED_LINE_SEARCH
85    DIVERGED_INNER              = SNES_DIVERGED_INNER
86    DIVERGED_LOCAL_MIN          = SNES_DIVERGED_LOCAL_MIN
87    DIVERGED_DTOL               = SNES_DIVERGED_DTOL
88    DIVERGED_TR_DELTA           = SNES_DIVERGED_TR_DELTA
89
90
91class SNESNewtonALCorrectionType(object):
92    """SNESNEWTONAL correction type.
93
94    See Also
95    --------
96    petsc.SNESNewtonALCorrectionType
97
98    """
99    EXACT = SNES_NEWTONAL_CORRECTION_EXACT
100    NORMAL = SNES_NEWTONAL_CORRECTION_NORMAL
101
102# --------------------------------------------------------------------
103
104
105cdef class SNES(Object):
106    """Nonlinear equations solver.
107
108    SNES is described in the `PETSc manual <petsc:manual/snes>`.
109
110    See Also
111    --------
112    petsc.SNES
113
114    """
115
116    Type = SNESType
117    NormSchedule = SNESNormSchedule
118    ConvergedReason = SNESConvergedReason
119    NewtonALCorrectionType = SNESNewtonALCorrectionType
120
121    # --- xxx ---
122
123    def __cinit__(self):
124        self.obj  = <PetscObject*> &self.snes
125        self.snes = NULL
126
127    # --- xxx ---
128
129    def view(self, Viewer viewer=None) -> None:
130        """View the solver.
131
132        Collective.
133
134        Parameters
135        ----------
136        viewer
137            A `Viewer` instance or `None` for the default viewer.
138
139        See Also
140        --------
141        petsc.SNESView
142
143        """
144        cdef PetscViewer cviewer = NULL
145        if viewer is not None: cviewer = viewer.vwr
146        CHKERR(SNESView(self.snes, cviewer))
147
148    def destroy(self) -> Self:
149        """Destroy the solver.
150
151        Collective.
152
153        See Also
154        --------
155        petsc.SNESDestroy
156
157        """
158        CHKERR(SNESDestroy(&self.snes))
159        return self
160
161    def create(self, comm: Comm | None = None) -> Self:
162        """Create a SNES solver.
163
164        Collective.
165
166        Parameters
167        ----------
168        comm
169            MPI communicator, defaults to `Sys.getDefaultComm`.
170
171        See Also
172        --------
173        Sys.getDefaultComm, petsc.SNESCreate
174
175        """
176        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_DEFAULT)
177        cdef PetscSNES newsnes = NULL
178        CHKERR(SNESCreate(ccomm, &newsnes))
179        CHKERR(PetscCLEAR(self.obj)); self.snes = newsnes
180        return self
181
182    def setType(self, snes_type: Type | str) -> None:
183        """Set the type of the solver.
184
185        Logically collective.
186
187        Parameters
188        ----------
189        snes_type
190            The type of the solver.
191
192        See Also
193        --------
194        getType, petsc.SNESSetType
195
196        """
197        cdef PetscSNESType cval = NULL
198        snes_type = str2bytes(snes_type, &cval)
199        CHKERR(SNESSetType(self.snes, cval))
200
201    def getType(self) -> str:
202        """Return the type of the solver.
203
204        Not collective.
205
206        See Also
207        --------
208        setType, petsc.SNESGetType
209
210        """
211        cdef PetscSNESType cval = NULL
212        CHKERR(SNESGetType(self.snes, &cval))
213        return bytes2str(cval)
214
215    def setOptionsPrefix(self, prefix: str | None) -> None:
216        """Set the prefix used for searching for options in the database.
217
218        Logically collective.
219
220        See Also
221        --------
222        petsc_options, petsc.SNESSetOptionsPrefix
223
224        """
225        cdef const char *cval = NULL
226        prefix = str2bytes(prefix, &cval)
227        CHKERR(SNESSetOptionsPrefix(self.snes, cval))
228
229    def getOptionsPrefix(self) -> str:
230        """Return the prefix used for searching for options in the database.
231
232        Not collective.
233
234        See Also
235        --------
236        petsc_options, setOptionsPrefix, petsc.SNESGetOptionsPrefix
237
238        """
239        cdef const char *cval = NULL
240        CHKERR(SNESGetOptionsPrefix(self.snes, &cval))
241        return bytes2str(cval)
242
243    def appendOptionsPrefix(self, prefix: str | None) -> None:
244        """Append to the prefix used for searching for options in the database.
245
246        Logically collective.
247
248        See Also
249        --------
250        petsc_options, setOptionsPrefix, petsc.SNESAppendOptionsPrefix
251
252        """
253        cdef const char *cval = NULL
254        prefix = str2bytes(prefix, &cval)
255        CHKERR(SNESAppendOptionsPrefix(self.snes, cval))
256
257    def setFromOptions(self) -> None:
258        """Configure the solver from the options database.
259
260        Collective.
261
262        See Also
263        --------
264        petsc_options, petsc.SNESSetFromOptions
265
266        """
267        CHKERR(SNESSetFromOptions(self.snes))
268
269    # --- application context ---
270
271    def setApplicationContext(self, appctx: Any) -> None:
272        """Set the application context."""
273        self.set_attr('__appctx__', appctx)
274        if appctx is not None:
275            registerAppCtx(<void*>appctx)
276            CHKERR(SNESSetApplicationContext(self.snes, <void*>appctx))
277        else:
278            CHKERR(SNESSetApplicationContext(self.snes, NULL))
279
280    def getApplicationContext(self) -> Any:
281        """Return the application context."""
282        cdef void *ctx = NULL
283        appctx = self.get_attr('__appctx__')
284        if appctx is None:
285            CHKERR(SNESGetApplicationContext(self.snes, &ctx))
286            appctx = toAppCtx(ctx)
287        return appctx
288
289    # backward compatibility
290    setAppCtx = setApplicationContext
291    getAppCtx = getApplicationContext
292
293    # --- discretization space ---
294
295    def getDM(self) -> DM:
296        """Return the `DM` associated with the solver.
297
298        Not collective.
299
300        See Also
301        --------
302        setDM, petsc.SNESGetDM
303
304        """
305        cdef PetscDM newdm = NULL
306        CHKERR(SNESGetDM(self.snes, &newdm))
307        cdef DM dm = subtype_DM(newdm)()
308        dm.dm = newdm
309        CHKERR(PetscINCREF(dm.obj))
310        return dm
311
312    def setDM(self, DM dm) -> None:
313        """Associate a `DM` with the solver.
314
315        Not collective.
316
317        See Also
318        --------
319        getDM, petsc.SNESSetDM
320
321        """
322        CHKERR(SNESSetDM(self.snes, dm.dm))
323
324    # --- TR ---
325
326    def setTRTolerances(self, delta_min: float | None = None, delta_max: float | None = None, delta_0: float | None = None) -> None:
327        """Set the tolerance parameters used for the trust region.
328
329        Logically collective.
330
331        Parameters
332        ----------
333        delta_min
334            The minimum allowed trust region size. Defaults to `CURRENT`.
335        delta_max
336            The maximum allowed trust region size. Defaults to `CURRENT`.
337        delta_0
338            The initial trust region size. Defaults to `CURRENT`.
339
340        See Also
341        --------
342        getTRTolerances, setTRUpdateParameters
343        petsc.SNESNewtonTRSetTolerances
344
345        """
346        cdef PetscReal cdmin, cdmax, cd0
347        cdmin = cdmax = cd0 = PETSC_CURRENT
348        if delta_min is not None: cdmin  = asReal(delta_min)
349        if delta_max is not None: cdmax  = asReal(delta_max)
350        if delta_0   is not None: cd0    = asReal(delta_0)
351        CHKERR(SNESNewtonTRSetTolerances(self.snes, cdmin, cdmax, cd0))
352
353    def getTRTolerances(self) -> tuple[float, float, float]:
354        """Return the tolerance parameters used for the trust region.
355
356        Not collective.
357
358        Returns
359        -------
360        delta_min : float
361            The minimum allowed trust region size.
362        delta_max : float
363            The maximum allowed trust region size.
364        delta_0 : float
365            The initial trust region size.
366
367        See Also
368        --------
369        setTRTolerances, getTRUpdateParameters
370        petsc.SNESNewtonTRGetTolerances
371
372        """
373        cdef PetscReal cdmin = 0, cdmax = 0, cd0 = 0
374        CHKERR(SNESNewtonTRGetTolerances(self.snes, &cdmin, &cdmax, &cd0))
375        return (toReal(cdmin), toReal(cdmax), toReal(cd0))
376
377    def setTRUpdateParameters(self,
378                              eta1: float | None = None, eta2: float | None = None, eta3: float | None = None,
379                              t1: float | None = None, t2: float | None = None) -> None:
380        """Set the update parameters used for the trust region.
381
382        Logically collective.
383
384        Parameters
385        ----------
386        eta1
387            The step acceptance tolerance. Defaults to `CURRENT`.
388        eta2
389            The shrinking tolerance. Defaults to `CURRENT`.
390        eta3
391            The enlarging tolerance. Defaults to `CURRENT`.
392        t1
393            The shrinking factor. Defaults to `CURRENT`.
394        t2
395            The enlarging factor. Defaults to `CURRENT`.
396
397        See Also
398        --------
399        setTRTolerances, getTRUpdateParameters
400        petsc.SNESNewtonTRSetUpdateParameters
401
402        """
403        cdef PetscReal ceta1, ceta2, ceta3, ct1, ct2
404        ceta1 = ceta2 = ceta3 = ct1 = ct2 = PETSC_CURRENT
405        if eta1 is not None: ceta1 = asReal(eta1)
406        if eta2 is not None: ceta2 = asReal(eta2)
407        if eta3 is not None: ceta3 = asReal(eta3)
408        if t1 is not None: ct1 = asReal(t1)
409        if t2 is not None: ct2 = asReal(t2)
410        CHKERR(SNESNewtonTRSetUpdateParameters(self.snes, ceta1, ceta2, ceta3, ct1, ct2))
411
412    def getTRUpdateParameters(self) -> tuple[float, float, float, float, float]:
413        """Return the update parameters used for the trust region.
414
415        Not collective.
416
417        Returns
418        -------
419        eta1 : float
420            The step acceptance tolerance.
421        eta2 : float
422            The shrinking tolerance.
423        eta3 : float
424            The enlarging tolerance.
425        t1 : float
426            The shrinking factor.
427        t2 : float
428            The enlarging factor.
429
430        See Also
431        --------
432        setTRUpdateParameters, getTRTolerances
433        petsc.SNESNewtonTRGetUpdateParameters
434
435        """
436        cdef PetscReal ceta1 = 0, ceta2 = 0, ceta3 = 0, ct1 = 0, ct2 = 0
437        CHKERR(SNESNewtonTRGetUpdateParameters(self.snes, &ceta1, &ceta2, &ceta3, &ct1, &ct2))
438        return (toReal(ceta1), toReal(ceta2), toReal(ceta3), toReal(ct1), toReal(ct2))
439
440    # --- FAS ---
441
442    def setFASInterpolation(self, level: int, Mat mat) -> None:
443        """Set the `Mat` to be used to apply the interpolation from level-1 to level.
444
445        Collective.
446
447        See Also
448        --------
449        getFASInterpolation, setFASRestriction, setFASInjection
450        petsc.SNESFASSetInterpolation, petsc.SNESFAS
451
452        """
453        cdef PetscInt clevel = asInt(level)
454        CHKERR(SNESFASSetInterpolation(self.snes, clevel, mat.mat))
455
456    def getFASInterpolation(self, level: int) -> Mat:
457        """Return the `Mat` used to apply the interpolation from level-1 to level.
458
459        Not collective.
460
461        See Also
462        --------
463        setFASInterpolation, petsc.SNESFASGetInterpolation, petsc.SNESFAS
464
465        """
466        cdef PetscInt clevel = asInt(level)
467        cdef Mat mat = Mat()
468        CHKERR(SNESFASGetInterpolation(self.snes, clevel, &mat.mat))
469        CHKERR(PetscINCREF(mat.obj))
470        return mat
471
472    def setFASRestriction(self, level: int, Mat mat) -> None:
473        """Set the `Mat` to be used to apply the restriction from level-1 to level.
474
475        Collective.
476
477        See Also
478        --------
479        setFASRScale, getFASRestriction, setFASInterpolation, setFASInjection
480        petsc.SNESFASSetRestriction, petsc.SNESFAS
481
482        """
483        cdef PetscInt clevel = asInt(level)
484        CHKERR(SNESFASSetRestriction(self.snes, clevel, mat.mat))
485
486    def getFASRestriction(self, level: int) -> Mat:
487        """Return the `Mat` used to apply the restriction from level-1 to level.
488
489        Not collective.
490
491        See Also
492        --------
493        setFASRestriction, petsc.SNESFASGetRestriction, petsc.SNESFAS
494
495        """
496        cdef PetscInt clevel = asInt(level)
497        cdef Mat mat = Mat()
498        CHKERR(SNESFASGetRestriction(self.snes, clevel, &mat.mat))
499        CHKERR(PetscINCREF(mat.obj))
500        return mat
501
502    def setFASInjection(self, level: int, Mat mat) -> None:
503        """Set the `Mat` to be used to apply the injection from level-1 to level.
504
505        Collective.
506
507        See Also
508        --------
509        getFASInjection, setFASInterpolation, setFASRestriction
510        petsc.SNESFASSetInjection, petsc.SNESFAS
511
512        """
513        cdef PetscInt clevel = asInt(level)
514        CHKERR(SNESFASSetInjection(self.snes, clevel, mat.mat))
515
516    def getFASInjection(self, level: int) -> Mat:
517        """Return the `Mat` used to apply the injection from level-1 to level.
518
519        Not collective.
520
521        See Also
522        --------
523        setFASInjection, petsc.SNESFASGetInjection, petsc.SNESFAS
524
525        """
526        cdef PetscInt clevel = asInt(level)
527        cdef Mat mat = Mat()
528        CHKERR(SNESFASGetInjection(self.snes, clevel, &mat.mat))
529        CHKERR(PetscINCREF(mat.obj))
530        return mat
531
532    def setFASRScale(self, level: int, Vec vec) -> None:
533        """Set the scaling factor of the restriction operator from level to level-1.
534
535        Collective.
536
537        See Also
538        --------
539        setFASRestriction, petsc.SNESFASSetRScale, petsc.SNESFAS
540
541        """
542        cdef PetscInt clevel = asInt(level)
543        CHKERR(SNESFASSetRScale(self.snes, clevel, vec.vec))
544
545    def setFASLevels(self, levels: int, comms: Sequence[Comm] | None = None) -> None:
546        """Set the number of levels to use with FAS.
547
548        Collective.
549
550        Parameters
551        ----------
552        levels
553            The number of levels
554        comms
555            An optional sequence of communicators of length `levels`,
556            or `None` for the default communicator `Sys.getDefaultComm`.
557
558        See Also
559        --------
560        getFASLevels, petsc.SNESFASSetLevels, petsc.SNESFAS
561
562        """
563        cdef PetscInt clevels = asInt(levels)
564        cdef MPI_Comm *ccomms = NULL
565        cdef Py_ssize_t i = 0
566        if comms is not None:
567            if clevels != <PetscInt>len(comms):
568                raise ValueError("Must provide as many communicators as levels")
569            CHKERR(PetscMalloc(sizeof(MPI_Comm)*<size_t>clevels, &ccomms))
570            try:
571                for i, comm in enumerate(comms):
572                    ccomms[i] = def_Comm(comm, MPI_COMM_NULL)
573                CHKERR(SNESFASSetLevels(self.snes, clevels, ccomms))
574            finally:
575                CHKERR(PetscFree(ccomms))
576        else:
577            CHKERR(SNESFASSetLevels(self.snes, clevels, ccomms))
578
579    def getFASLevels(self) -> int:
580        """Return the number of levels used.
581
582        Not collective.
583
584        See Also
585        --------
586        setFASLevels, petsc.SNESFASGetLevels, petsc.SNESFAS
587
588        """
589        cdef PetscInt levels = 0
590        CHKERR(SNESFASGetLevels(self.snes, &levels))
591        return toInt(levels)
592
593    def getFASCycleSNES(self, level: int) -> SNES:
594        """Return the `SNES` corresponding to a particular level of the FAS hierarchy.
595
596        Not collective.
597
598        See Also
599        --------
600        setFASLevels, getFASCoarseSolve, getFASSmoother
601        petsc.SNESFASGetCycleSNES, petsc.SNESFAS
602
603        """
604        cdef PetscInt clevel = asInt(level)
605        cdef SNES lsnes = SNES()
606        CHKERR(SNESFASGetCycleSNES(self.snes, clevel, &lsnes.snes))
607        CHKERR(PetscINCREF(lsnes.obj))
608        return lsnes
609
610    def getFASCoarseSolve(self) -> SNES:
611        """Return the `SNES` used at the coarsest level of the FAS hierarchy.
612
613        Not collective.
614
615        See Also
616        --------
617        getFASSmoother, petsc.SNESFASGetCoarseSolve, petsc.SNESFAS
618
619        """
620        cdef SNES smooth = SNES()
621        CHKERR(SNESFASGetCoarseSolve(self.snes, &smooth.snes))
622        CHKERR(PetscINCREF(smooth.obj))
623        return smooth
624
625    def getFASSmoother(self, level: int) -> SNES:
626        """Return the smoother used at a given level of the FAS hierarchy.
627
628        Not collective.
629
630        See Also
631        --------
632        setFASLevels, getFASCoarseSolve, getFASSmootherDown, getFASSmootherUp
633        petsc.SNESFASGetSmoother, petsc.SNESFAS
634
635        """
636        cdef PetscInt clevel = asInt(level)
637        cdef SNES smooth = SNES()
638        CHKERR(SNESFASGetSmoother(self.snes, clevel, &smooth.snes))
639        CHKERR(PetscINCREF(smooth.obj))
640        return smooth
641
642    def getFASSmootherDown(self, level: int) -> SNES:
643        """Return the downsmoother used at a given level of the FAS hierarchy.
644
645        Not collective.
646
647        See Also
648        --------
649        setFASLevels, getFASCoarseSolve, getFASSmoother, getFASSmootherUp
650        petsc.SNESFASGetSmootherDown, petsc.SNESFAS
651
652        """
653        cdef PetscInt clevel = asInt(level)
654        cdef SNES smooth = SNES()
655        CHKERR(SNESFASGetSmootherDown(self.snes, clevel, &smooth.snes))
656        CHKERR(PetscINCREF(smooth.obj))
657        return smooth
658
659    def getFASSmootherUp(self, level: int) -> SNES:
660        """Return the upsmoother used at a given level of the FAS hierarchy.
661
662        Not collective.
663
664        See Also
665        --------
666        setFASLevels, getFASCoarseSolve, getFASSmoother, getFASSmootherDown
667        petsc.SNESFASGetSmootherUp, petsc.SNESFAS
668
669        """
670        cdef PetscInt clevel = asInt(level)
671        cdef SNES smooth = SNES()
672        CHKERR(SNESFASGetSmootherUp(self.snes, clevel, &smooth.snes))
673        CHKERR(PetscINCREF(smooth.obj))
674        return smooth
675
676    # --- nonlinear preconditioner ---
677
678    def getNPC(self) -> SNES:
679        """Return the nonlinear preconditioner associated with the solver.
680
681        Not collective.
682
683        See Also
684        --------
685        setNPC, hasNPC, setNPCSide, getNPCSide, petsc.SNESGetNPC
686
687        """
688        cdef SNES snes = SNES()
689        CHKERR(SNESGetNPC(self.snes, &snes.snes))
690        CHKERR(PetscINCREF(snes.obj))
691        return snes
692
693    def hasNPC(self) -> bool:
694        """Return a boolean indicating whether the solver has a nonlinear preconditioner.
695
696        Not collective.
697
698        See Also
699        --------
700        setNPC, getNPC, setNPCSide, getNPCSide, petsc.SNESHasNPC
701
702        """
703        cdef PetscBool flag = PETSC_FALSE
704        CHKERR(SNESHasNPC(self.snes, &flag))
705        return toBool(flag)
706
707    def setNPC(self, SNES snes) -> None:
708        """Set the nonlinear preconditioner.
709
710        Logically collective.
711
712        See Also
713        --------
714        getNPC, hasNPC, setNPCSide, getNPCSide, petsc.SNESSetNPC
715
716        """
717        CHKERR(SNESSetNPC(self.snes, snes.snes))
718
719    def setNPCSide(self, side: PC.Side) -> None:
720        """Set the nonlinear preconditioning side.
721
722        Collective.
723
724        See Also
725        --------
726        setNPC, getNPC, hasNPC, getNPCSide, petsc.SNESSetNPCSide
727
728        """
729        CHKERR(SNESSetNPCSide(self.snes, side))
730
731    def getNPCSide(self) -> PC.Side:
732        """Return the nonlinear preconditioning side.
733
734        Not collective.
735
736        See Also
737        --------
738        setNPC, getNPC, hasNPC, setNPCSide, petsc.SNESGetNPCSide
739
740        """
741        cdef PetscPCSide side = PC_RIGHT
742        CHKERR(SNESGetNPCSide(self.snes, &side))
743        return side
744
745    # --- user Function/Jacobian routines ---
746
747    def setLineSearchPreCheck(self, precheck: SNESLSPreFunction | None,
748                              args: tuple[Any, ...] | None = None,
749                              kargs: dict[str, Any] | None = None) -> None:
750        """Set the callback that will be called before applying the linesearch.
751
752        Logically collective.
753
754        Parameters
755        ----------
756        precheck
757            The callback.
758        args
759            Positional arguments for the callback.
760        kargs
761            Keyword arguments for the callback.
762
763        See Also
764        --------
765        petsc.SNESLineSearchSetPreCheck
766
767        """
768        cdef PetscSNESLineSearch snesls = NULL
769        CHKERR(SNESGetLineSearch(self.snes, &snesls))
770        if precheck is not None:
771            if args  is None: args  = ()
772            if kargs is None: kargs = {}
773            context = (precheck, args, kargs)
774            self.set_attr('__precheck__', context)
775            # FIXME callback
776            CHKERR(SNESLineSearchSetPreCheck(snesls, SNES_PreCheck, <void*> context))
777        else:
778            self.set_attr('__precheck__', None)
779            CHKERR(SNESLineSearchSetPreCheck(snesls, NULL, NULL))
780
781    def setInitialGuess(self, initialguess: SNESGuessFunction | None,
782                        args: tuple[Any, ...] | None = None,
783                        kargs: dict[str, Any] | None = None) -> None:
784        """Set the callback to compute the initial guess.
785
786        Logically collective.
787
788        Parameters
789        ----------
790        initialguess
791            The callback.
792        args
793            Positional arguments for the callback.
794        kargs
795            Keyword arguments for the callback.
796
797        See Also
798        --------
799        getInitialGuess, petsc.SNESSetComputeInitialGuess
800
801        """
802        if initialguess is not None:
803            if args  is None: args  = ()
804            if kargs is None: kargs = {}
805            context = (initialguess, args, kargs)
806            self.set_attr('__initialguess__', context)
807            CHKERR(SNESSetComputeInitialGuess(self.snes, SNES_InitialGuess, <void*>context))
808        else:
809            self.set_attr('__initialguess__', None)
810            CHKERR(SNESSetComputeInitialGuess(self.snes, NULL, NULL))
811
812    def getInitialGuess(self) -> SNESGuessFunction:
813        """Return the callback to compute the initial guess.
814
815        Not collective.
816
817        See Also
818        --------
819        setInitialGuess
820
821        """
822        return self.get_attr('__initialguess__')
823
824    def setFunction(self, function: SNESFunction | None,
825                    Vec f=None,
826                    args: tuple[Any, ...] | None = None,
827                    kargs: dict[str, Any] | None = None) -> None:
828        """Set the callback to compute the nonlinear function.
829
830        Logically collective.
831
832        Parameters
833        ----------
834        function
835            The callback.
836        f
837            An optional vector to store the result.
838        args
839            Positional arguments for the callback.
840        kargs
841            Keyword arguments for the callback.
842
843        See Also
844        --------
845        getFunction, petsc.SNESSetFunction
846
847        """
848        cdef PetscVec fvec=NULL
849        if f is not None: fvec = f.vec
850        if function is not None:
851            if args  is None: args  = ()
852            if kargs is None: kargs = {}
853            context = (function, args, kargs)
854            self.set_attr('__function__', context)
855            CHKERR(SNESSetFunction(self.snes, fvec, SNES_Function, <void*>context))
856        else:
857            CHKERR(SNESSetFunction(self.snes, fvec, NULL, NULL))
858
859    def getFunction(self) -> SNESFunction:
860        """Return the callback to compute the nonlinear function.
861
862        Not collective.
863
864        See Also
865        --------
866        setFunction, petsc.SNESGetFunction
867
868        """
869        cdef PetscErrorCode(*fun)(PetscSNES, PetscVec, PetscVec, void*) except PETSC_ERR_PYTHON nogil
870        cdef Vec f = Vec()
871        cdef void* ctx = NULL
872        fun = SNES_Function
873        CHKERR(SNESGetFunction(self.snes, &f.vec, &fun, &ctx))
874        CHKERR(PetscINCREF(f.obj))
875        cdef object function = self.get_attr('__function__')
876        cdef object context
877
878        if function is not None:
879            return (f, function)
880
881        if ctx != NULL and <void*>SNES_Function == <void*>fun:
882            context = <object>ctx
883            if context is not None:
884                assert type(context) is tuple
885                return (f, context)
886
887        return (f, None)
888
889    def setUpdate(self, update: SNESUpdateFunction | None,
890                  args: tuple[Any, ...] | None = None,
891                  kargs: dict[str, Any] | None = None) -> None:
892        """Set the callback to compute update at the beginning of each step.
893
894        Logically collective.
895
896        Parameters
897        ----------
898        update
899            The callback.
900        args
901            Positional arguments for the callback.
902        kargs
903            Keyword arguments for the callback.
904
905        See Also
906        --------
907        getUpdate, petsc.SNESSetUpdate
908
909        """
910        if update is not None:
911            if args  is None: args  = ()
912            if kargs is None: kargs = {}
913            context = (update, args, kargs)
914            self.set_attr('__update__', context)
915            CHKERR(SNESSetUpdate(self.snes, SNES_Update))
916        else:
917            self.set_attr('__update__', None)
918            CHKERR(SNESSetUpdate(self.snes, NULL))
919
920    def getUpdate(self) -> SNESUpdateFunction:
921        """Return the callback to compute the update at the beginning of each step.
922
923        Not collective.
924
925        See Also
926        --------
927        setUpdate
928
929        """
930        return self.get_attr('__update__')
931
932    def setJacobian(self,
933                    jacobian: SNESJacobianFunction | None,
934                    Mat J=None, Mat P=None,
935                    args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
936        """Set the callback to compute the Jacobian.
937
938        Logically collective.
939
940        Parameters
941        ----------
942        jacobian
943            The Jacobian callback.
944        J
945            The matrix to store the Jacobian.
946        P
947            The matrix to construct the preconditioner.
948        args
949            Positional arguments for the callback.
950        kargs
951            Keyword arguments for the callback.
952
953        See Also
954        --------
955        getJacobian, petsc.SNESSetJacobian
956
957        """
958        cdef PetscMat Jmat=NULL
959        if J is not None: Jmat = J.mat
960        cdef PetscMat Pmat=Jmat
961        if P is not None: Pmat = P.mat
962        if jacobian is not None:
963            if args  is None: args  = ()
964            if kargs is None: kargs = {}
965            context = (jacobian, args, kargs)
966            self.set_attr('__jacobian__', context)
967            CHKERR(SNESSetJacobian(self.snes, Jmat, Pmat, SNES_Jacobian, <void*>context))
968        else:
969            CHKERR(SNESSetJacobian(self.snes, Jmat, Pmat, NULL, NULL))
970
971    def getJacobian(self) -> tuple[Mat, Mat, SNESJacobianFunction]:
972        """Return the matrices used to compute the Jacobian and the callback tuple.
973
974        Not collective.
975
976        Returns
977        -------
978        J : Mat
979            The matrix to store the Jacobian.
980        P : Mat
981            The matrix to construct the preconditioner.
982        callback : SNESJacobianFunction
983            callback, positional and keyword arguments.
984
985        See Also
986        --------
987        setJacobian, petsc.SNESGetJacobian
988
989        """
990        cdef Mat J = Mat()
991        cdef Mat P = Mat()
992        CHKERR(SNESGetJacobian(self.snes, &J.mat, &P.mat, NULL, NULL))
993        CHKERR(PetscINCREF(J.obj))
994        CHKERR(PetscINCREF(P.obj))
995        cdef object jacobian = self.get_attr('__jacobian__')
996        return (J, P, jacobian)
997
998    def setObjective(self,
999                     objective: SNESObjFunction | None,
1000                     args: tuple[Any, ...] | None = None,
1001                     kargs: dict[str, Any] | None = None) -> None:
1002        """Set the callback to compute the objective function.
1003
1004        Logically collective.
1005
1006        Parameters
1007        ----------
1008        objective
1009            The Jacobian callback.
1010        args
1011            Positional arguments for the callback.
1012        kargs
1013            Keyword arguments for the callback.
1014
1015        See Also
1016        --------
1017        getObjective, petsc.SNESSetObjective
1018
1019        """
1020        if objective is not None:
1021            if args  is None: args  = ()
1022            if kargs is None: kargs = {}
1023            context = (objective, args, kargs)
1024            self.set_attr('__objective__', context)
1025            CHKERR(SNESSetObjective(self.snes, SNES_Objective, <void*>context))
1026        else:
1027            CHKERR(SNESSetObjective(self.snes, NULL, NULL))
1028
1029    def getObjective(self) -> SNESObjFunction:
1030        """Return the objective callback tuple.
1031
1032        Not collective.
1033
1034        See Also
1035        --------
1036        setObjective
1037
1038        """
1039        CHKERR(SNESGetObjective(self.snes, NULL, NULL))
1040        cdef object objective = self.get_attr('__objective__')
1041        return objective
1042
1043    def computeFunction(self, Vec x, Vec f) -> None:
1044        """Compute the function.
1045
1046        Collective.
1047
1048        Parameters
1049        ----------
1050        x
1051            The input state vector.
1052        f
1053            The output vector.
1054
1055        See Also
1056        --------
1057        setFunction, petsc.SNESComputeFunction
1058
1059        """
1060        CHKERR(SNESComputeFunction(self.snes, x.vec, f.vec))
1061
1062    def computeJacobian(self, Vec x, Mat J, Mat P=None) -> None:
1063        """Compute the Jacobian.
1064
1065        Collective.
1066
1067        Parameters
1068        ----------
1069        x
1070            The input state vector.
1071        J
1072            The output Jacobian matrix.
1073        P
1074            The output Jacobian matrix used to construct the preconditioner.
1075
1076        See Also
1077        --------
1078        setJacobian, petsc.SNESComputeJacobian
1079
1080        """
1081        cdef PetscMat jmat = J.mat, pmat = J.mat
1082        if P is not None: pmat = P.mat
1083        CHKERR(SNESComputeJacobian(self.snes, x.vec, jmat, pmat))
1084
1085    def computeObjective(self, Vec x) -> float:
1086        """Compute the value of the objective function.
1087
1088        Collective.
1089
1090        Parameters
1091        ----------
1092        x
1093            The input state vector.
1094
1095        See Also
1096        --------
1097        setObjective, petsc.SNESComputeObjective
1098
1099        """
1100        cdef PetscReal o = 0
1101        CHKERR(SNESComputeObjective(self.snes, x.vec, &o))
1102        return toReal(o)
1103
1104    def setNGS(self,
1105               ngs: SNESNGSFunction | None,
1106               args: tuple[Any, ...] | None = None,
1107               kargs: dict[str, Any] | None = None) -> None:
1108        """Set the callback to compute nonlinear Gauss-Seidel.
1109
1110        Logically collective.
1111
1112        Parameters
1113        ----------
1114        ngs
1115            The nonlinear Gauss-Seidel callback.
1116        args
1117            Positional arguments for the callback.
1118        kargs
1119            Keyword arguments for the callback.
1120
1121        See Also
1122        --------
1123        getNGS, computeNGS, petsc.SNESSetNGS
1124
1125        """
1126        if ngs is not None:
1127            if args  is None: args  = ()
1128            if kargs is None: kargs = {}
1129            context = (ngs, args, kargs)
1130            self.set_attr('__ngs__', context)
1131            CHKERR(SNESSetNGS(self.snes, SNES_NGS, <void*>context))
1132        else:
1133            CHKERR(SNESSetNGS(self.snes, NULL, NULL))
1134
1135    def getNGS(self) -> SNESNGSFunction:
1136        """Return the nonlinear Gauss-Seidel callback tuple.
1137
1138        Not collective.
1139
1140        See Also
1141        --------
1142        setNGS, computeNGS
1143
1144        """
1145        CHKERR(SNESGetNGS(self.snes, NULL, NULL))
1146        cdef object ngs = self.get_attr('__ngs__')
1147        return ngs
1148
1149    def computeNGS(self, Vec x, Vec b=None) -> None:
1150        """Compute a nonlinear Gauss-Seidel step.
1151
1152        Collective.
1153
1154        Parameters
1155        ----------
1156        x
1157            The input/output state vector.
1158        b
1159            The input right-hand side vector.
1160
1161        See Also
1162        --------
1163        setNGS, getNGS, petsc.SNESComputeNGS
1164
1165        """
1166        cdef PetscVec bvec = NULL
1167        if b is not None: bvec = b.vec
1168        CHKERR(SNESComputeNGS(self.snes, bvec, x.vec))
1169
1170    # --- tolerances and convergence ---
1171
1172    def setTolerances(self, rtol: float | None = None, atol: float | None = None, stol: float | None = None, max_it: int | None = None) -> None:
1173        """Set the tolerance parameters used in the solver convergence tests.
1174
1175        Logically collective.
1176
1177        Parameters
1178        ----------
1179        rtol
1180            The relative norm of the residual, or `DETERMINE` to
1181            use the value when the object's type was set. Defaults to
1182            `CURRENT`.
1183        atol
1184            The absolute norm of the residual, or `DETERMINE` to
1185            use the value when the object's type was set. Defaults to
1186            `CURRENT`.
1187        stol
1188            The absolute norm of the step, or `DETERMINE` to use
1189            the value when the object's type was set. Defaults to
1190            `CURRENT`.
1191        max_it
1192            The maximum allowed number of iterations, or `DETERMINE`
1193            to use the value when the object's type was set. Defaults
1194            to `CURRENT`. Use `UNLIMITED` to have no maximum.
1195
1196        See Also
1197        --------
1198        getTolerances, petsc.SNESSetTolerances
1199
1200        """
1201        cdef PetscReal crtol, catol, cstol
1202        crtol = catol = cstol = PETSC_CURRENT
1203        cdef PetscInt cmaxit = PETSC_CURRENT
1204        if rtol   is not None: crtol  = asReal(rtol)
1205        if atol   is not None: catol  = asReal(atol)
1206        if stol   is not None: cstol  = asReal(stol)
1207        if max_it is not None: cmaxit = asInt(max_it)
1208        CHKERR(SNESSetTolerances(self.snes, catol, crtol, cstol,
1209                                 cmaxit, PETSC_CURRENT))
1210
1211    def getTolerances(self) -> tuple[float, float, float, int]:
1212        """Return the tolerance parameters used in the solver convergence tests.
1213
1214        Not collective.
1215
1216        Returns
1217        -------
1218        rtol : float
1219            The relative norm of the residual.
1220        atol : float
1221            The absolute norm of the residual.
1222        stol : float
1223            The absolute norm of the step.
1224        max_it : int
1225            The maximum allowed number of iterations.
1226
1227        See Also
1228        --------
1229        setTolerances, petsc.SNESGetTolerances
1230
1231        """
1232        cdef PetscReal crtol=0, catol=0, cstol=0
1233        cdef PetscInt cmaxit=0
1234        CHKERR(SNESGetTolerances(self.snes, &catol, &crtol, &cstol,
1235                                 &cmaxit, NULL))
1236        return (toReal(crtol), toReal(catol), toReal(cstol), toInt(cmaxit))
1237
1238    def setDivergenceTolerance(self, dtol: float) -> None:
1239        """Set the divergence tolerance parameter used in the convergence tests.
1240
1241        Logically collective.
1242
1243        Parameters
1244        ----------
1245        dtol
1246            The divergence tolerance.
1247
1248        See Also
1249        --------
1250        getDivergenceTolerance, setTolerances, petsc.SNESSetDivergenceTolerance
1251
1252        """
1253        CHKERR(SNESSetDivergenceTolerance(self.snes, asReal(dtol)))
1254
1255    def getDivergenceTolerance(self) -> float:
1256        """Get the divergence tolerance parameter used in the convergence tests.
1257
1258        Not collective.
1259
1260        See Also
1261        --------
1262        setDivergenceTolerance, getTolerances, petsc.SNESGetDivergenceTolerance
1263
1264        """
1265        cdef PetscReal cdtol=0
1266        CHKERR(SNESGetDivergenceTolerance(self.snes, &cdtol))
1267        return toReal(cdtol)
1268
1269    def setNormSchedule(self, normsched: NormSchedule) -> None:
1270        """Set the norm schedule.
1271
1272        Collective.
1273
1274        See Also
1275        --------
1276        getNormSchedule, petsc.SNESSetNormSchedule
1277
1278        """
1279        CHKERR(SNESSetNormSchedule(self.snes, normsched))
1280
1281    def getNormSchedule(self) -> NormSchedule:
1282        """Return the norm schedule.
1283
1284        Not collective.
1285
1286        See Also
1287        --------
1288        setNormSchedule, petsc.SNESGetNormSchedule
1289
1290        """
1291        cdef PetscSNESNormSchedule normsched = SNES_NORM_NONE
1292        CHKERR(SNESGetNormSchedule(self.snes, &normsched))
1293        return normsched
1294
1295    def setConvergenceTest(self, converged: SNESConvergedFunction | Literal["skip", "default"],
1296                           args: tuple[Any, ...] | None = None,
1297                           kargs: dict[str, Any] | None = None) -> None:
1298        """Set the callback to use as convergence test.
1299
1300        Logically collective.
1301
1302        Parameters
1303        ----------
1304        converged
1305            The convergence testing callback.
1306        args
1307            Positional arguments for the callback.
1308        kargs
1309            Keyword arguments for the callback.
1310
1311        See Also
1312        --------
1313        getConvergenceTest, callConvergenceTest, petsc.SNESSetConvergenceTest
1314
1315        """
1316        if converged == "skip":
1317            self.set_attr('__converged__', None)
1318            CHKERR(SNESSetConvergenceTest(self.snes, SNESConvergedSkip, NULL, NULL))
1319        elif converged is None or converged == "default":
1320            self.set_attr('__converged__', None)
1321            CHKERR(SNESSetConvergenceTest(self.snes, SNESConvergedDefault, NULL, NULL))
1322        else:
1323            assert callable(converged)
1324            if args  is None: args  = ()
1325            if kargs is None: kargs = {}
1326            context = (converged, args, kargs)
1327            self.set_attr('__converged__', context)
1328            CHKERR(SNESSetConvergenceTest(self.snes, SNES_Converged, <void*>context, NULL))
1329
1330    def getConvergenceTest(self) -> SNESConvergedFunction:
1331        """Return the callback to used as convergence test.
1332
1333        Not collective.
1334
1335        See Also
1336        --------
1337        setConvergenceTest, callConvergenceTest
1338
1339        """
1340        return self.get_attr('__converged__')
1341
1342    def callConvergenceTest(self, its: int, xnorm: float, ynorm: float, fnorm: float) -> ConvergedReason:
1343        """Compute the convergence test.
1344
1345        Collective.
1346
1347        Parameters
1348        ----------
1349        its
1350            Iteration number.
1351        xnorm
1352            Solution norm.
1353        ynorm
1354            Update norm.
1355        fnorm
1356            Function norm.
1357
1358        See Also
1359        --------
1360        setConvergenceTest, getConvergenceTest
1361
1362        """
1363        cdef PetscInt  ival  = asInt(its)
1364        cdef PetscReal rval1 = asReal(xnorm)
1365        cdef PetscReal rval2 = asReal(ynorm)
1366        cdef PetscReal rval3 = asReal(fnorm)
1367        cdef PetscSNESConvergedReason reason = SNES_CONVERGED_ITERATING
1368        CHKERR(SNESConvergenceTestCall(self.snes, ival,
1369                                       rval1, rval2, rval3, &reason))
1370        return reason
1371
1372    def converged(self, its: int, xnorm: float, ynorm: float, fnorm: float) -> None:
1373        """Compute the convergence test and update the solver converged reason.
1374
1375        Collective.
1376
1377        Parameters
1378        ----------
1379        its
1380            Iteration number.
1381        xnorm
1382            Solution norm.
1383        ynorm
1384            Update norm.
1385        fnorm
1386            Function norm.
1387
1388        See Also
1389        --------
1390        setConvergenceTest, getConvergenceTest, petsc.SNESConverged
1391
1392        """
1393        cdef PetscInt  ival  = asInt(its)
1394        cdef PetscReal rval1 = asReal(xnorm)
1395        cdef PetscReal rval2 = asReal(ynorm)
1396        cdef PetscReal rval3 = asReal(fnorm)
1397        CHKERR(SNESConverged(self.snes, ival, rval1, rval2, rval3))
1398
1399    def setConvergenceHistory(self, length=None, reset=False) -> None:
1400        """Set the convergence history.
1401
1402        Logically collective.
1403
1404        See Also
1405        --------
1406        petsc.SNESSetConvergenceHistory
1407
1408        """
1409        cdef PetscReal *rdata = NULL
1410        cdef PetscInt  *idata = NULL
1411        cdef PetscInt   size = 1000
1412        cdef PetscBool flag = PETSC_FALSE
1413        # FIXME
1414        if   length is True:     pass
1415        elif length is not None: size = asInt(length)
1416        if size < 0: size = 1000
1417        if reset: flag = PETSC_TRUE
1418        cdef object rhist = oarray_r(empty_r(size), NULL, &rdata)
1419        cdef object ihist = oarray_i(empty_i(size), NULL, &idata)
1420        self.set_attr('__history__', (rhist, ihist))
1421        CHKERR(SNESSetConvergenceHistory(self.snes, rdata, idata, size, flag))
1422
1423    def getConvergenceHistory(self) -> tuple[ArrayReal, ArrayInt]:
1424        """Return the convergence history.
1425
1426        Not collective.
1427
1428        See Also
1429        --------
1430        petsc.SNESGetConvergenceHistory
1431
1432        """
1433        cdef PetscReal *rdata = NULL
1434        cdef PetscInt  *idata = NULL
1435        cdef PetscInt   size = 0
1436        CHKERR(SNESGetConvergenceHistory(self.snes, &rdata, &idata, &size))
1437        cdef object rhist = array_r(size, rdata)
1438        cdef object ihist = array_i(size, idata)
1439        return (rhist, ihist)
1440
1441    def logConvergenceHistory(self, norm: float, linear_its: int = 0) -> None:
1442        """Log residual norm and linear iterations."""
1443        cdef PetscReal rval = asReal(norm)
1444        cdef PetscInt  ival = asInt(linear_its)
1445        CHKERR(SNESLogConvergenceHistory(self.snes, rval, ival))
1446
1447    def setResetCounters(self, reset: bool = True) -> None:
1448        """Set the flag to reset the counters.
1449
1450        Collective.
1451
1452        See Also
1453        --------
1454        petsc.SNESSetCountersReset
1455
1456        """
1457        cdef PetscBool flag = reset
1458        CHKERR(SNESSetCountersReset(self.snes, flag))
1459
1460    # --- monitoring ---
1461
1462    def setMonitor(self,
1463                   monitor: SNESMonitorFunction | None,
1464                   args: tuple[Any, ...] | None = None,
1465                   kargs: dict[str, Any] | None = None) -> None:
1466        """Set the callback used to monitor solver convergence.
1467
1468        Logically collective.
1469
1470        Parameters
1471        ----------
1472        monitor
1473            The callback.
1474        args
1475            Positional arguments for the callback.
1476        kargs
1477            Keyword arguments for the callback.
1478
1479        See Also
1480        --------
1481        getMonitor, petsc.SNESMonitorSet
1482
1483        """
1484        if monitor is None: return
1485        cdef object monitorlist = self.get_attr('__monitor__')
1486        if monitorlist is None:
1487            monitorlist = []
1488            self.set_attr('__monitor__', monitorlist)
1489            CHKERR(SNESMonitorSet(self.snes, SNES_Monitor, NULL, NULL))
1490        if args  is None: args  = ()
1491        if kargs is None: kargs = {}
1492        context = (monitor, args, kargs)
1493        monitorlist.append(context)
1494
1495    def getMonitor(self) -> list[tuple[SNESMonitorFunction, tuple[Any, ...], dict[str, Any]]]:
1496        """Return the callback used to monitor solver convergence.
1497
1498        Not collective.
1499
1500        See Also
1501        --------
1502        setMonitor
1503
1504        """
1505        return self.get_attr('__monitor__')
1506
1507    def monitorCancel(self) -> None:
1508        """Cancel all the monitors of the solver.
1509
1510        Logically collective.
1511
1512        See Also
1513        --------
1514        setMonitor, petsc.SNESMonitorCancel
1515
1516        """
1517        CHKERR(SNESMonitorCancel(self.snes))
1518        self.set_attr('__monitor__', None)
1519
1520    cancelMonitor = monitorCancel
1521
1522    def monitor(self, its, rnorm) -> None:
1523        """Monitor the solver.
1524
1525        Collective.
1526
1527        Parameters
1528        ----------
1529        its
1530            Current number of iterations.
1531        rnorm
1532            Current value of the residual norm.
1533
1534        See Also
1535        --------
1536        setMonitor, petsc.SNESMonitor
1537
1538        """
1539        cdef PetscInt  ival = asInt(its)
1540        cdef PetscReal rval = asReal(rnorm)
1541        CHKERR(SNESMonitor(self.snes, ival, rval))
1542
1543    # --- more tolerances ---
1544
1545    def setMaxFunctionEvaluations(self, max_funcs: int) -> None:
1546        """Set the maximum allowed number of function evaluations.
1547
1548        Collective.
1549
1550        See Also
1551        --------
1552        getMaxFunctionEvaluations, petsc.SNESSetTolerances
1553
1554        """
1555        cdef PetscReal r = PETSC_CURRENT
1556        cdef PetscInt  i = PETSC_CURRENT
1557        cdef PetscInt ival = asInt(max_funcs)
1558        CHKERR(SNESSetTolerances(self.snes, r, r, r, i, ival))
1559
1560    def getMaxFunctionEvaluations(self) -> int:
1561        """Return the maximum allowed number of function evaluations.
1562
1563        Not collective.
1564
1565        See Also
1566        --------
1567        setMaxFunctionEvaluations, petsc.SNESSetTolerances
1568
1569        """
1570        cdef PetscReal *r = NULL
1571        cdef PetscInt  *i = NULL
1572        cdef PetscInt ival = 0
1573        CHKERR(SNESGetTolerances(self.snes, r, r, r, i, &ival))
1574        return toInt(ival)
1575
1576    def getFunctionEvaluations(self) -> int:
1577        """Return the current number of function evaluations.
1578
1579        Not collective.
1580
1581        See Also
1582        --------
1583        setMaxFunctionEvaluations, petsc.SNESGetNumberFunctionEvals
1584
1585        """
1586        cdef PetscInt ival = 0
1587        CHKERR(SNESGetNumberFunctionEvals(self.snes, &ival))
1588        return toInt(ival)
1589
1590    def setMaxStepFailures(self, max_fails: int) -> None:
1591        """Set the maximum allowed number of step failures.
1592
1593        Collective.
1594
1595        See Also
1596        --------
1597        getMaxStepFailures, petsc.SNESSetMaxNonlinearStepFailures
1598
1599        """
1600        cdef PetscInt ival = asInt(max_fails)
1601        CHKERR(SNESSetMaxNonlinearStepFailures(self.snes, ival))
1602
1603    def getMaxStepFailures(self) -> int:
1604        """Return the maximum allowed number of step failures.
1605
1606        Not collective.
1607
1608        See Also
1609        --------
1610        setMaxStepFailures, petsc.SNESGetMaxNonlinearStepFailures
1611
1612        """
1613        cdef PetscInt ival = 0
1614        CHKERR(SNESGetMaxNonlinearStepFailures(self.snes, &ival))
1615        return toInt(ival)
1616
1617    def getStepFailures(self) -> int:
1618        """Return the current number of step failures.
1619
1620        Not collective.
1621
1622        See Also
1623        --------
1624        getMaxStepFailures, petsc.SNESGetNonlinearStepFailures
1625
1626        """
1627        cdef PetscInt ival = 0
1628        CHKERR(SNESGetNonlinearStepFailures(self.snes, &ival))
1629        return toInt(ival)
1630
1631    def setMaxKSPFailures(self, max_fails: int) -> None:
1632        """Set the maximum allowed number of linear solve failures.
1633
1634        Collective.
1635
1636        See Also
1637        --------
1638        getMaxKSPFailures, petsc.SNESSetMaxLinearSolveFailures
1639
1640        """
1641        cdef PetscInt ival = asInt(max_fails)
1642        CHKERR(SNESSetMaxLinearSolveFailures(self.snes, ival))
1643
1644    def getMaxKSPFailures(self) -> int:
1645        """Return the maximum allowed number of linear solve failures.
1646
1647        Not collective.
1648
1649        See Also
1650        --------
1651        setMaxKSPFailures, petsc.SNESGetMaxLinearSolveFailures
1652
1653        """
1654        cdef PetscInt ival = 0
1655        CHKERR(SNESGetMaxLinearSolveFailures(self.snes, &ival))
1656        return toInt(ival)
1657
1658    def getKSPFailures(self) -> int:
1659        """Return the current number of linear solve failures.
1660
1661        Not collective.
1662
1663        See Also
1664        --------
1665        getMaxKSPFailures, petsc.SNESGetLinearSolveFailures
1666
1667        """
1668        cdef PetscInt ival = 0
1669        CHKERR(SNESGetLinearSolveFailures(self.snes, &ival))
1670        return toInt(ival)
1671
1672    setMaxNonlinearStepFailures = setMaxStepFailures
1673    getMaxNonlinearStepFailures = getMaxStepFailures
1674    getNonlinearStepFailures    = getStepFailures
1675    setMaxLinearSolveFailures   = setMaxKSPFailures
1676    getMaxLinearSolveFailures   = getMaxKSPFailures
1677    getLinearSolveFailures      = getKSPFailures
1678
1679    # --- solving ---
1680
1681    def setUp(self) -> None:
1682        """Set up the internal data structures for using the solver.
1683
1684        Collective.
1685
1686        See Also
1687        --------
1688        petsc.SNESSetUp
1689
1690        """
1691        CHKERR(SNESSetUp(self.snes))
1692
1693    def setUpMatrices(self) -> None:
1694        """Ensures that matrices are available for Newton-like methods.
1695
1696        Collective.
1697
1698        This is only of use to implementers of custom SNES types.
1699
1700        See Also
1701        --------
1702        setUp, petsc.SNESSetUpMatrices
1703
1704        """
1705        CHKERR(SNESSetUpMatrices(self.snes))
1706
1707    def reset(self) -> None:
1708        """Reset the solver.
1709
1710        Collective.
1711
1712        See Also
1713        --------
1714        petsc.SNESReset
1715
1716        """
1717        CHKERR(SNESReset(self.snes))
1718
1719    def solve(self, Vec b = None, Vec x = None) -> None:
1720        """Solve the nonlinear equations.
1721
1722        Collective.
1723
1724        Parameters
1725        ----------
1726        b
1727            The affine right-hand side or `None` to use zero.
1728        x
1729            The starting vector or `None` to use the vector stored internally.
1730
1731        See Also
1732        --------
1733        setSolution, getSolution, petsc.SNESSolve
1734
1735        """
1736        cdef PetscVec rhs = NULL
1737        cdef PetscVec sol = NULL
1738        if b is not None: rhs = b.vec
1739        if x is not None: sol = x.vec
1740        CHKERR(SNESSolve(self.snes, rhs, sol))
1741
1742    def setConvergedReason(self, reason: ConvergedReason) -> None:
1743        """Set the termination flag.
1744
1745        Collective.
1746
1747        See Also
1748        --------
1749        getConvergedReason, petsc.SNESSetConvergedReason
1750
1751        """
1752        cdef PetscSNESConvergedReason eval = reason
1753        CHKERR(SNESSetConvergedReason(self.snes, eval))
1754
1755    def getConvergedReason(self) -> ConvergedReason:
1756        """Return the termination flag.
1757
1758        Not collective.
1759
1760        See Also
1761        --------
1762        setConvergedReason, petsc.SNESGetConvergedReason
1763
1764        """
1765        cdef PetscSNESConvergedReason reason = SNES_CONVERGED_ITERATING
1766        CHKERR(SNESGetConvergedReason(self.snes, &reason))
1767        return reason
1768
1769    def setErrorIfNotConverged(self, flag: bool) -> None:
1770        """Immediately generate an error if the solver has not converged.
1771
1772        Collective.
1773
1774        See Also
1775        --------
1776        getErrorIfNotConverged, petsc.SNESSetErrorIfNotConverged
1777
1778        """
1779        cdef PetscBool ernc = asBool(flag)
1780        CHKERR(SNESSetErrorIfNotConverged(self.snes, ernc))
1781
1782    def getErrorIfNotConverged(self) -> bool:
1783        """Return the flag indicating error on divergence.
1784
1785        Not collective.
1786
1787        See Also
1788        --------
1789        setErrorIfNotConverged, petsc.SNESGetErrorIfNotConverged
1790
1791        """
1792        cdef PetscBool flag = PETSC_FALSE
1793        CHKERR(SNESGetErrorIfNotConverged(self.snes, &flag))
1794        return toBool(flag)
1795
1796    def setIterationNumber(self, its: int) -> None:
1797        """Set the current iteration number.
1798
1799        Collective.
1800
1801        This is only of use to implementers of custom SNES types.
1802
1803        See Also
1804        --------
1805        getIterationNumber, petsc.SNESSetIterationNumber
1806
1807        """
1808        cdef PetscInt ival = asInt(its)
1809        CHKERR(SNESSetIterationNumber(self.snes, ival))
1810
1811    def getIterationNumber(self) -> int:
1812        """Return the current iteration number.
1813
1814        Not collective.
1815
1816        See Also
1817        --------
1818        setIterationNumber, petsc.SNESGetIterationNumber
1819
1820        """
1821        cdef PetscInt ival = 0
1822        CHKERR(SNESGetIterationNumber(self.snes, &ival))
1823        return toInt(ival)
1824
1825    def setForceIteration(self, force: bool) -> None:
1826        """Force solve to take at least one iteration.
1827
1828        Collective.
1829
1830        See Also
1831        --------
1832        petsc.SNESSetForceIteration
1833
1834        """
1835        cdef PetscBool bval = asBool(force)
1836        CHKERR(SNESSetForceIteration(self.snes, bval))
1837
1838    def setFunctionNorm(self, norm: float) -> None:
1839        """Set the function norm value.
1840
1841        Collective.
1842
1843        This is only of use to implementers of custom SNES types.
1844
1845        See Also
1846        --------
1847        getFunctionNorm, petsc.SNESSetFunctionNorm
1848
1849        """
1850        cdef PetscReal rval = asReal(norm)
1851        CHKERR(SNESSetFunctionNorm(self.snes, rval))
1852
1853    def getFunctionNorm(self) -> float:
1854        """Return the function norm.
1855
1856        Not collective.
1857
1858        See Also
1859        --------
1860        setFunctionNorm, petsc.SNESGetFunctionNorm
1861
1862        """
1863        cdef PetscReal rval = 0
1864        CHKERR(SNESGetFunctionNorm(self.snes, &rval))
1865        return toReal(rval)
1866
1867    def getLinearSolveIterations(self) -> int:
1868        """Return the total number of linear iterations.
1869
1870        Not collective.
1871
1872        See Also
1873        --------
1874        petsc.SNESGetLinearSolveIterations
1875
1876        """
1877        cdef PetscInt ival = 0
1878        CHKERR(SNESGetLinearSolveIterations(self.snes, &ival))
1879        return toInt(ival)
1880
1881    def getRhs(self) -> Vec:
1882        """Return the vector holding the right-hand side.
1883
1884        Not collective.
1885
1886        See Also
1887        --------
1888        petsc.SNESGetRhs
1889
1890        """
1891        cdef Vec vec = Vec()
1892        CHKERR(SNESGetRhs(self.snes, &vec.vec))
1893        CHKERR(PetscINCREF(vec.obj))
1894        return vec
1895
1896    def getSolution(self) -> Vec:
1897        """Return the vector holding the solution.
1898
1899        Not collective.
1900
1901        See Also
1902        --------
1903        setSolution, petsc.SNESGetSolution
1904
1905        """
1906        cdef Vec vec = Vec()
1907        CHKERR(SNESGetSolution(self.snes, &vec.vec))
1908        CHKERR(PetscINCREF(vec.obj))
1909        return vec
1910
1911    def setSolution(self, Vec vec) -> None:
1912        """Set the vector used to store the solution.
1913
1914        Collective.
1915
1916        See Also
1917        --------
1918        getSolution, petsc.SNESSetSolution
1919
1920        """
1921        CHKERR(SNESSetSolution(self.snes, vec.vec))
1922
1923    def getSolutionUpdate(self) -> Vec:
1924        """Return the vector holding the solution update.
1925
1926        Not collective.
1927
1928        See Also
1929        --------
1930        petsc.SNESGetSolutionUpdate
1931
1932        """
1933        cdef Vec vec = Vec()
1934        CHKERR(SNESGetSolutionUpdate(self.snes, &vec.vec))
1935        CHKERR(PetscINCREF(vec.obj))
1936        return vec
1937
1938    # --- linear solver ---
1939
1940    def setKSP(self, KSP ksp) -> None:
1941        """Set the linear solver that will be used by the nonlinear solver.
1942
1943        Logically collective.
1944
1945        See Also
1946        --------
1947        getKSP, petsc.SNESSetKSP
1948
1949        """
1950        CHKERR(SNESSetKSP(self.snes, ksp.ksp))
1951
1952    def getKSP(self) -> KSP:
1953        """Return the linear solver used by the nonlinear solver.
1954
1955        Not collective.
1956
1957        See Also
1958        --------
1959        setKSP, petsc.SNESGetKSP
1960
1961        """
1962        cdef KSP ksp = KSP()
1963        CHKERR(SNESGetKSP(self.snes, &ksp.ksp))
1964        CHKERR(PetscINCREF(ksp.obj))
1965        return ksp
1966
1967    def setUseEW(self, flag: bool = True, *targs: Any, **kargs: Any) -> None:
1968        """Tell the solver to use the Eisenstat-Walker trick.
1969
1970        Logically collective.
1971
1972        Parameters
1973        ----------
1974        flag
1975            Whether or not to use the Eisenstat-Walker trick.
1976        *targs
1977            Positional arguments for `setParamsEW`.
1978        **kargs
1979            Keyword arguments for `setParamsEW`.
1980
1981        See Also
1982        --------
1983        getUseEW, setParamsEW, petsc.SNESKSPSetUseEW
1984
1985        """
1986        cdef PetscBool bval = flag
1987        CHKERR(SNESKSPSetUseEW(self.snes, bval))
1988        if targs or kargs: self.setParamsEW(*targs, **kargs)
1989
1990    def getUseEW(self) -> bool:
1991        """Return the flag indicating if the solver uses the Eisenstat-Walker trick.
1992
1993        Not collective.
1994
1995        See Also
1996        --------
1997        setUseEW, setParamsEW, petsc.SNESKSPGetUseEW
1998
1999        """
2000        cdef PetscBool flag = PETSC_FALSE
2001        CHKERR(SNESKSPGetUseEW(self.snes, &flag))
2002        return toBool(flag)
2003
2004    def setParamsEW(self,
2005                    version: int | None = None,
2006                    rtol_0: float | None = None,
2007                    rtol_max: float | None = None,
2008                    gamma: float | None = None,
2009                    alpha: float | None = None,
2010                    alpha2: float | None = None,
2011                    threshold: float | None = None) -> None:
2012        """Set the parameters for the Eisenstat and Walker trick.
2013
2014        Logically collective.
2015
2016        Parameters
2017        ----------
2018        version
2019            The version of the algorithm. Defaults to `CURRENT`.
2020        rtol_0
2021            The initial relative residual norm. Defaults to `CURRENT`.
2022        rtol_max
2023            The maximum relative residual norm. Defaults to `CURRENT`.
2024        gamma
2025            Parameter. Defaults to `CURRENT`.
2026        alpha
2027            Parameter. Defaults to `CURRENT`.
2028        alpha2
2029            Parameter. Defaults to `CURRENT`.
2030        threshold
2031            Parameter. Defaults to `CURRENT`.
2032
2033        See Also
2034        --------
2035        setUseEW, getParamsEW, petsc.SNESKSPSetParametersEW
2036
2037        """
2038        cdef PetscInt  cversion   = PETSC_CURRENT
2039        cdef PetscReal crtol_0    = PETSC_CURRENT
2040        cdef PetscReal crtol_max  = PETSC_CURRENT
2041        cdef PetscReal cgamma     = PETSC_CURRENT
2042        cdef PetscReal calpha     = PETSC_CURRENT
2043        cdef PetscReal calpha2    = PETSC_CURRENT
2044        cdef PetscReal cthreshold = PETSC_CURRENT
2045        if version   is not None: cversion   = asInt(version)
2046        if rtol_0    is not None: crtol_0    = asReal(rtol_0)
2047        if rtol_max  is not None: crtol_max  = asReal(rtol_max)
2048        if gamma     is not None: cgamma     = asReal(gamma)
2049        if alpha     is not None: calpha     = asReal(alpha)
2050        if alpha2    is not None: calpha2    = asReal(alpha2)
2051        if threshold is not None: cthreshold = asReal(threshold)
2052        CHKERR(SNESKSPSetParametersEW(
2053            self.snes, cversion, crtol_0, crtol_max,
2054            cgamma, calpha, calpha2, cthreshold))
2055
2056    def getParamsEW(self) -> dict[str, int | float]:
2057        """Get the parameters of the Eisenstat and Walker trick.
2058
2059        Not collective.
2060
2061        See Also
2062        --------
2063        setUseEW, setParamsEW, petsc.SNESKSPGetParametersEW
2064
2065        """
2066        cdef PetscInt  version=0
2067        cdef PetscReal rtol_0=0, rtol_max=0
2068        cdef PetscReal gamma=0, alpha=0, alpha2=0
2069        cdef PetscReal threshold=0
2070        CHKERR(SNESKSPGetParametersEW(
2071            self.snes, &version, &rtol_0, &rtol_max,
2072            &gamma, &alpha, &alpha2, &threshold))
2073        return {'version'   : toInt(version),
2074                'rtol_0'    : toReal(rtol_0),
2075                'rtol_max'  : toReal(rtol_max),
2076                'gamma'     : toReal(gamma),
2077                'alpha'     : toReal(alpha),
2078                'alpha2'    : toReal(alpha2),
2079                'threshold' : toReal(threshold), }
2080
2081    def setUseKSP(self, flag=True) -> None:
2082        """Set the boolean flag indicating to use a linear solver.
2083
2084        Logically collective.
2085
2086        See Also
2087        --------
2088        getUseKSP
2089
2090        """
2091        cdef PetscBool bval = flag
2092        CHKERR(SNESSetUseKSP(self.snes, bval))
2093
2094    def getUseKSP(self) -> bool:
2095        """Return the flag indicating if the solver uses a linear solver.
2096
2097        Not collective.
2098
2099        See Also
2100        --------
2101        setUseKSP
2102
2103        """
2104        cdef PetscBool flag = PETSC_FALSE
2105        CHKERR(SNESGetUseKSP(self.snes, &flag))
2106        return toBool(flag)
2107
2108    # --- matrix-free / finite differences ---
2109
2110    def setUseMF(self, flag=True) -> None:
2111        """Set the boolean flag indicating to use matrix-free finite-differencing.
2112
2113        Logically collective.
2114
2115        See Also
2116        --------
2117        getUseMF
2118
2119        """
2120        cdef PetscBool bval = flag
2121        CHKERR(SNESSetUseMFFD(self.snes, bval))
2122
2123    def getUseMF(self) -> bool:
2124        """Return the flag indicating if the solver uses matrix-free finite-differencing.
2125
2126        Not collective.
2127
2128        See Also
2129        --------
2130        setUseMF
2131
2132        """
2133        cdef PetscBool flag = PETSC_FALSE
2134        CHKERR(SNESGetUseMFFD(self.snes, &flag))
2135        return toBool(flag)
2136
2137    def setUseFD(self, flag=True) -> None:
2138        """Set the boolean flag to use coloring finite-differencing for Jacobian assembly.
2139
2140        Logically collective.
2141
2142        See Also
2143        --------
2144        getUseFD
2145
2146        """
2147        cdef PetscBool bval = flag
2148        CHKERR(SNESSetUseFDColoring(self.snes, bval))
2149
2150    def getUseFD(self) -> bool:
2151        """Return ``true`` if the solver uses color finite-differencing for the Jacobian.
2152
2153        Not collective.
2154
2155        See Also
2156        --------
2157        setUseFD
2158
2159        """
2160        cdef PetscBool flag = PETSC_FALSE
2161        CHKERR(SNESGetUseFDColoring(self.snes, &flag))
2162        return toBool(flag)
2163
2164    # --- VI ---
2165
2166    def setVariableBounds(self, Vec xl, Vec xu) -> None:
2167        """Set the vector for the variable bounds.
2168
2169        Collective.
2170
2171        See Also
2172        --------
2173        getVariableBounds, petsc.SNESVISetVariableBounds
2174
2175        """
2176        CHKERR(SNESVISetVariableBounds(self.snes, xl.vec, xu.vec))
2177
2178    def getVariableBounds(self) -> tuple[Vec, Vec]:
2179        """Get the vectors for the variable bounds.
2180
2181        Collective.
2182
2183        See Also
2184        --------
2185        setVariableBounds, petsc.SNESVIGetVariableBounds
2186
2187        """
2188        cdef Vec xl = Vec(), xu = Vec()
2189        CHKERR(SNESVIGetVariableBounds(self.snes, &xl.vec, &xu.vec))
2190        CHKERR(PetscINCREF(xl.obj)); CHKERR(PetscINCREF(xu.obj))
2191        return (xl, xu)
2192
2193    def getVIInactiveSet(self) -> IS:
2194        """Return the index set for the inactive set.
2195
2196        Not collective.
2197
2198        See Also
2199        --------
2200        petsc.SNESVIGetInactiveSet
2201
2202        """
2203        cdef IS inact = IS()
2204        CHKERR(SNESVIGetInactiveSet(self.snes, &inact.iset))
2205        CHKERR(PetscINCREF(inact.obj))
2206        return inact
2207
2208    # --- Python ---
2209
2210    def createPython(self, context: Any = None, comm: Comm | None = None) -> Self:
2211        """Create a nonlinear solver of Python type.
2212
2213        Collective.
2214
2215        Parameters
2216        ----------
2217        context
2218            An instance of the Python class implementing the required methods.
2219        comm
2220            MPI communicator, defaults to `Sys.getDefaultComm`.
2221
2222        See Also
2223        --------
2224        petsc_python_snes, setType, setPythonContext, Type.PYTHON
2225
2226        """
2227        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_DEFAULT)
2228        cdef PetscSNES newsnes = NULL
2229        CHKERR(SNESCreate(ccomm, &newsnes))
2230        CHKERR(PetscCLEAR(self.obj)); self.snes = newsnes
2231        CHKERR(SNESSetType(self.snes, SNESPYTHON))
2232        CHKERR(SNESPythonSetContext(self.snes, <void*>context))
2233        return self
2234
2235    def setPythonContext(self, context: Any) -> None:
2236        """Set the instance of the class implementing the required Python methods.
2237
2238        Not collective.
2239
2240        See Also
2241        --------
2242        petsc_python_snes, getPythonContext
2243
2244        """
2245        CHKERR(SNESPythonSetContext(self.snes, <void*>context))
2246
2247    def getPythonContext(self) -> Any:
2248        """Return the instance of the class implementing the required Python methods.
2249
2250        Not collective.
2251
2252        See Also
2253        --------
2254        petsc_python_snes, setPythonContext
2255
2256        """
2257        cdef void *context = NULL
2258        CHKERR(SNESPythonGetContext(self.snes, &context))
2259        if context == NULL: return None
2260        else: return <object> context
2261
2262    def setPythonType(self, py_type: str) -> None:
2263        """Set the fully qualified Python name of the class to be used.
2264
2265        Collective.
2266
2267        See Also
2268        --------
2269        petsc_python_snes, setPythonContext, getPythonType
2270        petsc.SNESPythonSetType
2271
2272        """
2273        cdef const char *cval = NULL
2274        py_type = str2bytes(py_type, &cval)
2275        CHKERR(SNESPythonSetType(self.snes, cval))
2276
2277    def getPythonType(self) -> str:
2278        """Return the fully qualified Python name of the class used by the solver.
2279
2280        Not collective.
2281
2282        See Also
2283        --------
2284        petsc_python_snes, setPythonContext, setPythonType
2285        petsc.SNESPythonGetType
2286
2287        """
2288        cdef const char *cval = NULL
2289        CHKERR(SNESPythonGetType(self.snes, &cval))
2290        return bytes2str(cval)
2291
2292    # --- Composite ---
2293
2294    def getCompositeSNES(self, int n) -> SNES:
2295        """Return the n-th solver in the composite.
2296
2297        Not collective.
2298
2299        See Also
2300        --------
2301        getCompositeNumber, petsc.SNESCompositeGetSNES, petsc.SNESCOMPOSITE
2302
2303        """
2304        cdef PetscInt cn
2305        cdef SNES snes = SNES()
2306        cn = asInt(n)
2307        CHKERR(SNESCompositeGetSNES(self.snes, cn, &snes.snes))
2308        CHKERR(PetscINCREF(snes.obj))
2309        return snes
2310
2311    def getCompositeNumber(self) -> int:
2312        """Return the number of solvers in the composite.
2313
2314        Not collective.
2315
2316        See Also
2317        --------
2318        getCompositeSNES, petsc.SNESCompositeGetNumber, petsc.SNESCOMPOSITE
2319
2320        """
2321        cdef PetscInt cn = 0
2322        CHKERR(SNESCompositeGetNumber(self.snes, &cn))
2323        return toInt(cn)
2324
2325    # --- NASM ---
2326
2327    def getNASMSNES(self, n: int) -> SNES:
2328        """Return the n-th solver in NASM.
2329
2330        Not collective.
2331
2332        See Also
2333        --------
2334        getNASMNumber, petsc.SNESNASMGetSNES, petsc.SNESNASM
2335
2336        """
2337        cdef PetscInt cn = asInt(n)
2338        cdef SNES snes = SNES()
2339        CHKERR(SNESNASMGetSNES(self.snes, cn, &snes.snes))
2340        CHKERR(PetscINCREF(snes.obj))
2341        return snes
2342
2343    def getNASMNumber(self) -> int:
2344        """Return the number of solvers in NASM.
2345
2346        Not collective.
2347
2348        See Also
2349        --------
2350        getNASMSNES, petsc.SNESNASMGetNumber, petsc.SNESNASM
2351
2352        """
2353        cdef PetscInt cn = 0
2354        CHKERR(SNESNASMGetNumber(self.snes, &cn))
2355        return toInt(cn)
2356
2357    # --- Patch ---
2358
2359    def setPatchCellNumbering(self, Section sec) -> None:
2360        """Set cell patch numbering."""
2361        CHKERR(SNESPatchSetCellNumbering(self.snes, sec.sec))
2362
2363    def setPatchDiscretisationInfo(self, dms, bs,
2364                                   cellNodeMaps,
2365                                   subspaceOffsets,
2366                                   ghostBcNodes,
2367                                   globalBcNodes) -> None:
2368        """Set patch discretisation information."""
2369        cdef PetscInt numSubSpaces = 0
2370        cdef PetscInt numGhostBcs = 0, numGlobalBcs = 0
2371        cdef PetscInt *nodesPerCell = NULL
2372        cdef const PetscInt **ccellNodeMaps = NULL
2373        cdef PetscDM *cdms = NULL
2374        cdef PetscInt *cbs = NULL
2375        cdef PetscInt *csubspaceOffsets = NULL
2376        cdef PetscInt *cghostBcNodes = NULL
2377        cdef PetscInt *cglobalBcNodes = NULL
2378        cdef PetscInt i = 0
2379
2380        bs = iarray_i(bs, &numSubSpaces, &cbs)
2381        ghostBcNodes = iarray_i(ghostBcNodes, &numGhostBcs, &cghostBcNodes)
2382        globalBcNodes = iarray_i(globalBcNodes, &numGlobalBcs, &cglobalBcNodes)
2383        subspaceOffsets = iarray_i(subspaceOffsets, NULL, &csubspaceOffsets)
2384
2385        CHKERR(PetscMalloc(<size_t>numSubSpaces*sizeof(PetscInt), &nodesPerCell))
2386        CHKERR(PetscMalloc(<size_t>numSubSpaces*sizeof(PetscDM), &cdms))
2387        CHKERR(PetscMalloc(<size_t>numSubSpaces*sizeof(PetscInt*), &ccellNodeMaps))
2388        for i in range(numSubSpaces):
2389            cdms[i] = (<DM?>dms[i]).dm
2390            _, nodes = asarray(cellNodeMaps[i]).shape
2391            cellNodeMaps[i] = iarray_i(cellNodeMaps[i], NULL, <PetscInt**>&ccellNodeMaps[i])
2392            nodesPerCell[i] = asInt(nodes)
2393
2394        # TODO: refactor on the PETSc side to take ISes?
2395        CHKERR(SNESPatchSetDiscretisationInfo(self.snes, numSubSpaces,
2396                                              cdms, cbs, nodesPerCell,
2397                                              ccellNodeMaps, csubspaceOffsets,
2398                                              numGhostBcs, cghostBcNodes,
2399                                              numGlobalBcs, cglobalBcNodes))
2400        CHKERR(PetscFree(nodesPerCell))
2401        CHKERR(PetscFree(cdms))
2402        CHKERR(PetscFree(ccellNodeMaps))
2403
2404    def setPatchComputeOperator(self, operator, args=None, kargs=None) -> None:
2405        """Set patch compute operator."""
2406        if args is  None: args  = ()
2407        if kargs is None: kargs = {}
2408        context = (operator, args, kargs)
2409        self.set_attr("__patch_compute_operator__", context)
2410        CHKERR(SNESPatchSetComputeOperator(self.snes, PCPatch_ComputeOperator, <void*>context))
2411
2412    def setPatchComputeFunction(self, function, args=None, kargs=None) -> None:
2413        """Set patch compute function."""
2414        if args is  None: args  = ()
2415        if kargs is None: kargs = {}
2416        context = (function, args, kargs)
2417        self.set_attr("__patch_compute_function__", context)
2418        CHKERR(SNESPatchSetComputeFunction(self.snes, PCPatch_ComputeFunction, <void*>context))
2419
2420    def setPatchConstructType(self, typ, operator=None, args=None, kargs=None) -> None:
2421        """Set patch construct type."""
2422        if args is  None: args  = ()
2423        if kargs is None: kargs = {}
2424
2425        if typ in {PC.PatchConstructType.PYTHON, PC.PatchConstructType.USER} and operator is None:
2426            raise ValueError("Must provide operator for USER or PYTHON type")
2427        if operator is not None:
2428            context = (operator, args, kargs)
2429        else:
2430            context = None
2431        self.set_attr("__patch_construction_operator__", context)
2432        CHKERR(SNESPatchSetConstructType(self.snes, typ, PCPatch_UserConstructOperator, <void*>context))
2433
2434    # --- application context ---
2435
2436    property appctx:
2437        """Application context."""
2438        def __get__(self) -> Any:
2439            return self.getAppCtx()
2440
2441        def __set__(self, value):
2442            self.setAppCtx(value)
2443
2444    # --- discretization space ---
2445
2446    property dm:
2447        """`DM`."""
2448        def __get__(self) -> DM:
2449            return self.getDM()
2450
2451        def __set__(self, value):
2452            self.setDM(value)
2453
2454    # --- nonlinear preconditioner ---
2455
2456    property npc:
2457        """Nonlinear preconditioner."""
2458        def __get__(self) -> SNES:
2459            return self.getNPC()
2460
2461        def __set__(self, value):
2462            self.setNPC(value)
2463
2464    # --- vectors ---
2465
2466    property vec_sol:
2467        """Solution vector."""
2468        def __get__(self) -> Vec:
2469            return self.getSolution()
2470
2471    property vec_upd:
2472        """Update vector."""
2473        def __get__(self) -> Vec:
2474            return self.getSolutionUpdate()
2475
2476    property vec_rhs:
2477        """Right-hand side vector."""
2478        def __get__(self) -> Vec:
2479            return self.getRhs()
2480
2481    # --- linear solver ---
2482
2483    property ksp:
2484        """Linear solver."""
2485        def __get__(self) -> KSP:
2486            return self.getKSP()
2487
2488        def __set__(self, value):
2489            self.setKSP(value)
2490
2491    property use_ksp:
2492        """Boolean indicating if the solver uses a linear solver."""
2493        def __get__(self) -> bool:
2494            return self.getUseKSP()
2495
2496        def __set__(self, value):
2497            self.setUseKSP(value)
2498
2499    property use_ew:
2500        """Use the Eisenstat-Walker trick."""
2501        def __get__(self) -> bool:
2502            return self.getUseEW()
2503
2504        def __set__(self, value):
2505            self.setUseEW(value)
2506
2507    # --- tolerances ---
2508
2509    property rtol:
2510        """Relative residual tolerance."""
2511        def __get__(self) -> float:
2512            return self.getTolerances()[0]
2513
2514        def __set__(self, value):
2515            self.setTolerances(rtol=value)
2516
2517    property atol:
2518        """Absolute residual tolerance."""
2519        def __get__(self) -> float:
2520            return self.getTolerances()[1]
2521
2522        def __set__(self, value):
2523            self.setTolerances(atol=value)
2524
2525    property stol:
2526        """Solution update tolerance."""
2527        def __get__(self) -> float:
2528            return self.getTolerances()[2]
2529
2530        def __set__(self, value):
2531            self.setTolerances(stol=value)
2532
2533    property max_it:
2534        """Maximum number of iterations."""
2535        def __get__(self) -> int:
2536            return self.getTolerances()[3]
2537
2538        def __set__(self, value):
2539            self.setTolerances(max_it=value)
2540
2541    # --- more tolerances ---
2542
2543    property max_funcs:
2544        """Maximum number of function evaluations."""
2545        def __get__(self) -> int:
2546            return self.getMaxFunctionEvaluations()
2547
2548        def __set__(self, value):
2549            self.setMaxFunctionEvaluations(value)
2550
2551    # --- iteration ---
2552
2553    property its:
2554        """Number of iterations."""
2555        def __get__(self) -> int:
2556            return self.getIterationNumber()
2557
2558        def __set__(self, value):
2559            self.setIterationNumber(value)
2560
2561    property norm:
2562        """Function norm."""
2563        def __get__(self) -> float:
2564            return self.getFunctionNorm()
2565
2566        def __set__(self, value):
2567            self.setFunctionNorm(value)
2568
2569    property history:
2570        """Convergence history."""
2571        def __get__(self) -> tuple[ArrayReal, ArrayInt]:
2572            return self.getConvergenceHistory()
2573
2574    # --- convergence ---
2575
2576    property reason:
2577        """Converged reason."""
2578        def __get__(self) -> ConvergedReason:
2579            return self.getConvergedReason()
2580
2581        def __set__(self, value):
2582            self.setConvergedReason(value)
2583
2584    property is_iterating:
2585        """Boolean indicating if the solver has not converged yet."""
2586        def __get__(self) -> bool:
2587            return self.reason == 0
2588
2589    property is_converged:
2590        """Boolean indicating if the solver has converged."""
2591        def __get__(self) -> bool:
2592            return self.reason > 0
2593
2594    property is_diverged:
2595        """Boolean indicating if the solver has failed."""
2596        def __get__(self) -> bool:
2597            return self.reason < 0
2598
2599    # --- matrix-free / finite differences ---
2600
2601    property use_mf:
2602        """Boolean indicating if the solver uses matrix-free finite-differencing."""
2603        def __get__(self) -> bool:
2604            return self.getUseMF()
2605
2606        def __set__(self, value):
2607            self.setUseMF(value)
2608
2609    property use_fd:
2610        """Boolean indicating if the solver uses coloring finite-differencing."""
2611        def __get__(self) -> bool:
2612            return self.getUseFD()
2613
2614        def __set__(self, value):
2615            self.setUseFD(value)
2616
2617    # --- linesearch ---
2618
2619    def getLineSearch(self) -> SNESLineSearch:
2620        """Return the linesearch object associated with this SNES.
2621
2622        Not collective.
2623
2624        See Also
2625        --------
2626        setLineSearch, linesearch, petsc.SNESGetLineSearch
2627
2628        """
2629        cdef SNESLineSearch linesearch = SNESLineSearch()
2630        CHKERR(SNESGetLineSearch(self.snes, &linesearch.snesls))
2631        CHKERR(PetscINCREF(linesearch.obj))
2632        return linesearch
2633
2634    def setLineSearch(self, SNESLineSearch linesearch) -> None:
2635        """Set the linesearch object to be used by this SNES.
2636
2637        Logically collective.
2638
2639        See Also
2640        --------
2641        getLineSearch, linesearch, petsc.SNESSetLineSearch
2642
2643        """
2644        CHKERR(SNESSetLineSearch(self.snes, linesearch.snesls))
2645
2646    property linesearch:
2647        """The linesearch object associated with this SNES."""
2648        def __get__(self) -> SNESLineSearch:
2649            return self.getLineSearch()
2650
2651        def __set__(self, value):
2652            self.setLineSearch(value)
2653
2654    # --- NewtonAL methods ---
2655
2656    def setNewtonALFunction(self, function: SNESFunction | None,
2657                            args: tuple[Any, ...] | None = None,
2658                            kargs: dict[str, Any] | None = None) -> None:
2659        """Set the callback to compute the tangent load vector for SNESNEWTONAL.
2660
2661        Logically collective.
2662
2663        Parameters
2664        ----------
2665        function
2666            The callback.
2667        args
2668            Positional arguments for the callback.
2669        kargs
2670            Keyword arguments for the callback.
2671
2672        See Also
2673        --------
2674        getNewtonALLoadParameter, petsc.SNESNewtonALSetFunction
2675
2676        """
2677        if function is not None:
2678            if args  is None: args  = ()
2679            if kargs is None: kargs = {}
2680            context = (function, args, kargs)
2681            self.set_attr('__newtonal_function__', context)
2682            CHKERR(SNESNewtonALSetFunction(self.snes, SNES_NewtonALFunction, <void*>context))
2683        else:
2684            self.set_attr('__newtonal_function__', None)
2685            CHKERR(SNESNewtonALSetFunction(self.snes, NULL, NULL))
2686
2687    def getNewtonALLoadParameter(self) -> float:
2688        """Return the load parameter for SNESNEWTONAL.
2689
2690        Not collective.
2691
2692        See Also
2693        --------
2694        setNewtonALFunction, setNewtonALCorrectionType
2695        petsc.SNESNewtonALGetLoadParameter
2696
2697        """
2698        cdef PetscReal load = 0
2699        CHKERR(SNESNewtonALGetLoadParameter(self.snes, &load))
2700        return toReal(load)
2701
2702    def setNewtonALCorrectionType(self, corrtype: NewtonALCorrectionType) -> None:
2703        """Set the correction type for SNESNEWTONAL.
2704
2705        Logically collective.
2706
2707        See Also
2708        --------
2709        getNewtonALLoadParameter, petsc.SNESNewtonALSetCorrectionType
2710
2711        """
2712        CHKERR(SNESNewtonALSetCorrectionType(self.snes, corrtype))
2713
2714# --------------------------------------------------------------------
2715
2716
2717class SNESLineSearchType(object):
2718    """SNES linesearch type.
2719
2720    See Also
2721    --------
2722    petsc.SNESLineSearchType
2723
2724    """
2725    BT         = S_(SNESLINESEARCHBT)
2726    NLEQERR    = S_(SNESLINESEARCHNLEQERR)
2727    BASIC      = S_(SNESLINESEARCHBASIC)
2728    NONE       = S_(SNESLINESEARCHNONE)
2729    SECANT     = S_(SNESLINESEARCHSECANT)
2730    CP         = S_(SNESLINESEARCHCP)
2731    SHELL      = S_(SNESLINESEARCHSHELL)
2732    NCGLINEAR  = S_(SNESLINESEARCHNCGLINEAR)
2733    BISECTION  = S_(SNESLINESEARCHBISECTION)
2734
2735
2736cdef class SNESLineSearch(Object):
2737    """Linesearch object used by SNES solvers.
2738
2739    See Also
2740    --------
2741    petsc.SNESLineSearch
2742
2743    """
2744
2745    Type = SNESLineSearchType
2746
2747    def __cinit__(self):
2748        self.obj  = <PetscObject*> &self.snesls
2749        self.snesls = NULL
2750
2751    def create(self, comm: Comm | None = None) -> Self:
2752        """Create a new linesearch object.
2753
2754        Collective.
2755
2756        Parameters
2757        ----------
2758        comm : Comm, optional
2759            MPI communicator, defaults to `Sys.getDefaultComm`.
2760
2761        See Also
2762        --------
2763        destroy, petsc.SNESLineSearchCreate
2764
2765        """
2766        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_DEFAULT)
2767        cdef PetscSNESLineSearch newls = NULL
2768        CHKERR(SNESLineSearchCreate(ccomm, &newls))
2769        CHKERR(PetscCLEAR(self.obj)); self.snesls = newls
2770        return self
2771
2772    def setFromOptions(self) -> None:
2773        """Configure the linesearch from the options database.
2774
2775        Collective.
2776
2777        See Also
2778        --------
2779        petsc.SNESLineSearchSetFromOptions
2780
2781        """
2782        CHKERR(SNESLineSearchSetFromOptions(self.snesls))
2783
2784    def view(self, Viewer viewer = None) -> None:
2785        """View the linesearch object.
2786
2787        Collective.
2788
2789        Parameters
2790        ----------
2791        viewer
2792            A `Viewer` instance or `None` for the default viewer.
2793
2794        See Also
2795        --------
2796        petsc.SNESLineSearchView
2797
2798        """
2799        cdef PetscViewer cviewer = NULL
2800        if viewer is not None: cviewer = viewer.vwr
2801        CHKERR(SNESLineSearchView(self.snesls, cviewer))
2802
2803    def getType(self) -> str:
2804        """Return the type of the linesearch.
2805
2806        Not collective.
2807
2808        See Also
2809        --------
2810        setType, petsc.SNESLineSearchGetType
2811
2812        """
2813        cdef const char *ctype = NULL
2814        CHKERR(SNESLineSearchGetType(self.snesls, &ctype))
2815        return bytes2str(ctype)
2816
2817    def setType(self, ls_type: Type | str) -> None:
2818        """Set the type of the linesearch.
2819
2820        Logically collective.
2821
2822        See Also
2823        --------
2824        getType, petsc.SNESLineSearchSetType
2825
2826        """
2827        cdef const char *ctype = NULL
2828        ls_type = str2bytes(ls_type, &ctype)
2829        CHKERR(SNESLineSearchSetType(self.snesls, ctype))
2830
2831    def getTolerances(self) -> tuple[float, float, float, float, float, int]:
2832        """Return the tolerance parameters used in the linesearch.
2833
2834        Not collective.
2835
2836        Returns
2837        -------
2838        minstep : float
2839            Minimum step length.
2840        maxstep : float, optional
2841            Maximum step length.
2842        rtol : float
2843            Relative tolerance.
2844        atol : float
2845            Absolute tolerance.
2846        ltol : float
2847            Lambda tolerance.
2848        max_its : int
2849            Maximum number of iterations for the linesearch.
2850
2851        See Also
2852        --------
2853        setTolerances, petsc.SNESLineSearchGetTolerances
2854        """
2855        cdef PetscReal rtol=0, atol=0, minstep=0, ltol=0, maxstep=0
2856        cdef PetscInt max_its=0
2857        CHKERR(SNESLineSearchGetTolerances(self.snesls, &minstep, &maxstep, &rtol, &atol, &ltol, &max_its))
2858        return (toReal(minstep), toReal(maxstep), toReal(rtol), toReal(atol), toReal(ltol), toInt(max_its))
2859
2860    def setTolerances(self, minstep: float | None = None, maxstep: float | None = None, rtol: float | None = None, atol: float | None = None, ltol: float | None = None, max_its: int | None = None) -> None:
2861        """Set the tolerance parameters used in the linesearch.
2862
2863        Logically collective.
2864
2865        Parameters
2866        ----------
2867        minstep : float
2868            Minimum step length.
2869        maxstep : float, optional
2870            Maximum step length.
2871        rtol : float, optional
2872            Relative tolerance.
2873        atol : float, optional
2874            Absolute tolerance.
2875        ltol : float, optional
2876            Lambda tolerance.
2877        max_its : int, optional
2878            Maximum number of iterations for the linesearch.
2879
2880        See Also
2881        --------
2882        getTolerances, petsc.SNESLineSearchSetTolerances
2883        """
2884        cdef PetscReal cminstep = PETSC_DEFAULT
2885        cdef PetscReal cmaxstep = PETSC_DEFAULT
2886        cdef PetscReal crtol = PETSC_DEFAULT
2887        cdef PetscReal catol = PETSC_DEFAULT
2888        cdef PetscReal cltol = PETSC_DEFAULT
2889        cdef PetscInt cmaxits = PETSC_DEFAULT
2890        if minstep is not None: cminstep = asReal(minstep)
2891        if maxstep is not None: cmaxstep = asReal(maxstep)
2892        if rtol is not None: crtol = asReal(rtol)
2893        if atol is not None: catol = asReal(atol)
2894        if ltol is not None: cltol = asReal(ltol)
2895        if max_its is not None: cmaxits = asInt(max_its)
2896        CHKERR(SNESLineSearchSetTolerances(self.snesls, cminstep, cmaxstep, crtol, catol, cltol, cmaxits))
2897
2898    def getOrder(self) -> int:
2899        """Return the order of the linesearch.
2900
2901        Not collective.
2902
2903        See Also
2904        --------
2905        setOrder, petsc.SNESLineSearchGetOrder
2906
2907        """
2908        cdef PetscInt order = 0
2909        CHKERR(SNESLineSearchGetOrder(self.snesls, &order))
2910        return toInt(order)
2911
2912    def setOrder(self, order: int) -> None:
2913        """Set the order of the linesearch.
2914
2915        Logically collective.
2916
2917        See Also
2918        --------
2919        getOrder, petsc.SNESLineSearchSetOrder
2920
2921        """
2922        cdef PetscInt iorder = asInt(order)
2923        CHKERR(SNESLineSearchSetOrder(self.snesls, iorder))
2924
2925    def destroy(self) -> Self:
2926        """Destroy the linesearch object.
2927
2928        Collective.
2929
2930        See Also
2931        --------
2932        petsc.SNESLineSearchDestroy
2933
2934        """
2935        CHKERR(SNESLineSearchDestroy(&self.snesls))
2936        return self
2937
2938del SNESType
2939del SNESNormSchedule
2940del SNESConvergedReason
2941del SNESNewtonALCorrectionType
2942del SNESLineSearchType
2943
2944# --------------------------------------------------------------------
2945