xref: /petsc/src/binding/petsc4py/src/petsc4py/PETSc/TAO.pyx (revision 124b60a56262a80503e09b9eaaec281e19388b1e) !
1# --------------------------------------------------------------------
2
3class TAOType:
4    """TAO solver type.
5
6    See Also
7    --------
8    petsc.TaoType
9
10    """
11    LMVM     = S_(TAOLMVM)
12    NLS      = S_(TAONLS)
13    NTR      = S_(TAONTR)
14    NTL      = S_(TAONTL)
15    CG       = S_(TAOCG)
16    TRON     = S_(TAOTRON)
17    OWLQN    = S_(TAOOWLQN)
18    BMRM     = S_(TAOBMRM)
19    BLMVM    = S_(TAOBLMVM)
20    BQNLS    = S_(TAOBQNLS)
21    BNCG     = S_(TAOBNCG)
22    BNLS     = S_(TAOBNLS)
23    BNTR     = S_(TAOBNTR)
24    BNTL     = S_(TAOBNTL)
25    BQNKLS   = S_(TAOBQNKLS)
26    BQNKTR   = S_(TAOBQNKTR)
27    BQNKTL   = S_(TAOBQNKTL)
28    BQPIP    = S_(TAOBQPIP)
29    GPCG     = S_(TAOGPCG)
30    NM       = S_(TAONM)
31    POUNDERS = S_(TAOPOUNDERS)
32    BRGN     = S_(TAOBRGN)
33    LCL      = S_(TAOLCL)
34    SSILS    = S_(TAOSSILS)
35    SSFLS    = S_(TAOSSFLS)
36    ASILS    = S_(TAOASILS)
37    ASFLS    = S_(TAOASFLS)
38    IPM      = S_(TAOIPM)
39    PDIPM    = S_(TAOPDIPM)
40    SHELL    = S_(TAOSHELL)
41    ADMM     = S_(TAOADMM)
42    ALMM     = S_(TAOALMM)
43    PYTHON   = S_(TAOPYTHON)
44
45
46class TAOConvergedReason:
47    """TAO solver termination reason.
48
49    See Also
50    --------
51    petsc.TaoConvergedReason
52
53    """
54    # iterating
55    CONTINUE_ITERATING    = TAO_CONTINUE_ITERATING    # iterating
56    CONVERGED_ITERATING   = TAO_CONTINUE_ITERATING    # iterating
57    ITERATING             = TAO_CONTINUE_ITERATING    # iterating
58    # converged
59    CONVERGED_GATOL       = TAO_CONVERGED_GATOL       # ||g(X)|| < gatol
60    CONVERGED_GRTOL       = TAO_CONVERGED_GRTOL       # ||g(X)||/f(X)  < grtol
61    CONVERGED_GTTOL       = TAO_CONVERGED_GTTOL       # ||g(X)||/||g(X0)|| < gttol
62    CONVERGED_STEPTOL     = TAO_CONVERGED_STEPTOL     # small step size
63    CONVERGED_MINF        = TAO_CONVERGED_MINF        # f(X) < F_min
64    CONVERGED_USER        = TAO_CONVERGED_USER        # user defined
65    # diverged
66    DIVERGED_MAXITS       = TAO_DIVERGED_MAXITS       #
67    DIVERGED_NAN          = TAO_DIVERGED_NAN          #
68    DIVERGED_MAXFCN       = TAO_DIVERGED_MAXFCN       #
69    DIVERGED_LS_FAILURE   = TAO_DIVERGED_LS_FAILURE   #
70    DIVERGED_TR_REDUCTION = TAO_DIVERGED_TR_REDUCTION #
71    DIVERGED_USER         = TAO_DIVERGED_USER         # user defined
72
73
74class TAOBNCGType:
75    """TAO Bound Constrained Conjugate Gradient (BNCG) Update Type."""
76    GD         = TAO_BNCG_GD
77    PCGD       = TAO_BNCG_PCGD
78    HS         = TAO_BNCG_HS
79    FR         = TAO_BNCG_FR
80    PRP        = TAO_BNCG_PRP
81    PRP_PLUS   = TAO_BNCG_PRP_PLUS
82    DY         = TAO_BNCG_DY
83    HZ         = TAO_BNCG_HZ
84    DK         = TAO_BNCG_DK
85    KD         = TAO_BNCG_KD
86    SSML_BFGS  = TAO_BNCG_SSML_BFGS
87    SSML_DFP   = TAO_BNCG_SSML_DFP
88    SSML_BRDN  = TAO_BNCG_SSML_BRDN
89
90
91class TAOALMMType:
92    """TAO Augmented Lagrangian Multiplier method (ALMM) Type."""
93    CLASSIC = TAO_ALMM_CLASSIC
94    PHR     = TAO_ALMM_PHR
95# --------------------------------------------------------------------
96
97
98cdef class TAO(Object):
99    """Optimization solver.
100
101    TAO is described in the `PETSc manual <petsc:manual/tao>`.
102
103    See Also
104    --------
105    petsc.Tao
106
107    """
108
109    Type = TAOType
110    ConvergedReason = TAOConvergedReason
111    BNCGType = TAOBNCGType
112    ALMMType = TAOALMMType
113    # FIXME backward compatibility
114    Reason = TAOConvergedReason
115
116    def __cinit__(self):
117        self.obj = <PetscObject*> &self.tao
118        self.tao = NULL
119
120    def view(self, Viewer viewer=None) -> None:
121        """View the solver.
122
123        Collective.
124
125        Parameters
126        ----------
127        viewer
128            A `Viewer` instance or `None` for the default viewer.
129
130        See Also
131        --------
132        petsc.TaoView
133
134        """
135        cdef PetscViewer vwr = NULL
136        if viewer is not None: vwr = viewer.vwr
137        CHKERR(TaoView(self.tao, vwr))
138
139    def destroy(self) -> Self:
140        """Destroy the solver.
141
142        Collective.
143
144        See Also
145        --------
146        petsc.TaoDestroy
147
148        """
149        CHKERR(TaoDestroy(&self.tao))
150        return self
151
152    def create(self, comm: Comm | None = None) -> Self:
153        """Create a TAO solver.
154
155        Collective.
156
157        Parameters
158        ----------
159        comm
160            MPI communicator, defaults to `Sys.getDefaultComm`.
161
162        See Also
163        --------
164        Sys.getDefaultComm, petsc.TaoCreate
165
166        """
167        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_DEFAULT)
168        cdef PetscTAO newtao = NULL
169        CHKERR(TaoCreate(ccomm, &newtao))
170        CHKERR(PetscCLEAR(self.obj)); self.tao = newtao
171        return self
172
173    def setType(self, tao_type: Type | str) -> None:
174        """Set the type of the solver.
175
176        Logically collective.
177
178        Parameters
179        ----------
180        tao_type
181            The type of the solver.
182
183        See Also
184        --------
185        getType, petsc.TaoSetType
186
187        """
188        cdef PetscTAOType ctype = NULL
189        tao_type = str2bytes(tao_type, &ctype)
190        CHKERR(TaoSetType(self.tao, ctype))
191
192    def getType(self) -> str:
193        """Return the type of the solver.
194
195        Not collective.
196
197        See Also
198        --------
199        setType, petsc.TaoGetType
200
201        """
202        cdef PetscTAOType ctype = NULL
203        CHKERR(TaoGetType(self.tao, &ctype))
204        return bytes2str(ctype)
205
206    def setOptionsPrefix(self, prefix: str | None) -> None:
207        """Set the prefix used for searching for options in the database.
208
209        Logically collective.
210
211        See Also
212        --------
213        petsc_options, petsc.TaoSetOptionsPrefix
214
215        """
216        cdef const char *cprefix = NULL
217        prefix = str2bytes(prefix, &cprefix)
218        CHKERR(TaoSetOptionsPrefix(self.tao, cprefix))
219
220    def appendOptionsPrefix(self, prefix: str | None) -> None:
221        """Append to the prefix used for searching for options in the database.
222
223        Logically collective.
224
225        See Also
226        --------
227        petsc_options, setOptionsPrefix, petsc.TaoAppendOptionsPrefix
228
229        """
230        cdef const char *cprefix = NULL
231        prefix = str2bytes(prefix, &cprefix)
232        CHKERR(TaoAppendOptionsPrefix(self.tao, cprefix))
233
234    def getOptionsPrefix(self) -> str:
235        """Return the prefix used for searching for options in the database.
236
237        Not collective.
238
239        See Also
240        --------
241        petsc_options, setOptionsPrefix, petsc.TaoGetOptionsPrefix
242
243        """
244        cdef const char *prefix = NULL
245        CHKERR(TaoGetOptionsPrefix(self.tao, &prefix))
246        return bytes2str(prefix)
247
248    def setFromOptions(self) -> None:
249        """Configure the solver from the options database.
250
251        Collective.
252
253        See Also
254        --------
255        petsc_options, petsc.TaoSetFromOptions
256
257        """
258        CHKERR(TaoSetFromOptions(self.tao))
259
260    def setUp(self) -> None:
261        """Set up the internal data structures for using the solver.
262
263        Collective.
264
265        See Also
266        --------
267        petsc.TaoSetUp
268
269        """
270        CHKERR(TaoSetUp(self.tao))
271
272    #
273
274    def setInitialTrustRegionRadius(self, radius: float) -> None:
275        """Set the initial trust region radius.
276
277        Collective.
278
279        See Also
280        --------
281        petsc.TaoSetInitialTrustRegionRadius
282
283        """
284        cdef PetscReal cradius = asReal(radius)
285        CHKERR(TaoSetInitialTrustRegionRadius(self.tao, cradius))
286
287    # --------------
288
289    def setAppCtx(self, appctx: Any) -> None:
290        """Set the application context."""
291        self.set_attr("__appctx__", appctx)
292
293    def getAppCtx(self) -> Any:
294        """Return the application context."""
295        return self.get_attr("__appctx__")
296
297    def setSolution(self, Vec x) -> None:
298        """Set the vector used to store the solution.
299
300        Collective.
301
302        See Also
303        --------
304        getSolution, petsc.TaoSetSolution
305
306        """
307        CHKERR(TaoSetSolution(self.tao, x.vec))
308
309    def setObjective(self, objective: TAOObjectiveFunction, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
310        """Set the objective function evaluation callback.
311
312        Logically collective.
313
314        Parameters
315        ----------
316        objective
317            The objective function callback.
318        args
319            Positional arguments for the callback.
320        kargs
321            Keyword arguments for the callback.
322
323        See Also
324        --------
325        setGradient, setObjectiveGradient, petsc.TaoSetObjective
326
327        """
328        if args is None: args = ()
329        if kargs is None: kargs = {}
330        context = (objective, args, kargs)
331        self.set_attr("__objective__", context)
332        CHKERR(TaoSetObjective(self.tao, TAO_Objective, <void*>context))
333
334    def setResidual(self, residual: TAOResidualFunction, Vec R, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
335        """Set the residual evaluation callback for least-squares applications.
336
337        Logically collective.
338
339        Parameters
340        ----------
341        residual
342            The residual callback.
343        R
344            The vector to store the residual.
345        args
346            Positional arguments for the callback.
347        kargs
348            Keyword arguments for the callback.
349
350        See Also
351        --------
352        setJacobianResidual, petsc.TaoSetResidualRoutine
353
354        """
355        if args is None: args = ()
356        if kargs is None: kargs = {}
357        context = (residual, args, kargs)
358        self.set_attr("__residual__", context)
359        CHKERR(TaoSetResidualRoutine(self.tao, R.vec, TAO_Residual, <void*>context))
360
361    def setJacobianResidual(self, jacobian: TAOJacobianResidualFunction, Mat J=None, Mat P=None, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
362        """Set the callback to compute the least-squares residual Jacobian.
363
364        Logically collective.
365
366        Parameters
367        ----------
368        jacobian
369            The Jacobian callback.
370        J
371            The matrix to store the Jacobian.
372        P
373            The matrix to construct the preconditioner.
374        args
375            Positional arguments for the callback.
376        kargs
377            Keyword arguments for the callback.
378
379        See Also
380        --------
381        setResidual, petsc.TaoSetJacobianResidualRoutine
382
383        """
384        cdef PetscMat Jmat = NULL
385        if J is not None: Jmat = J.mat
386        cdef PetscMat Pmat = Jmat
387        if P is not None: Pmat = P.mat
388        if args is None: args = ()
389        if kargs is None: kargs = {}
390        context = (jacobian, args, kargs)
391        self.set_attr("__jacobian_residual__", context)
392        CHKERR(TaoSetJacobianResidualRoutine(self.tao, Jmat, Pmat, TAO_JacobianResidual, <void*>context))
393
394    def setGradient(self, gradient: TAOGradientFunction, Vec g=None, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
395        """Set the gradient evaluation callback.
396
397        Logically collective.
398
399        Parameters
400        ----------
401        gradient
402            The gradient callback.
403        g
404            The vector to store the gradient.
405        args
406            Positional arguments for the callback.
407        kargs
408            Keyword arguments for the callback.
409
410        See Also
411        --------
412        setObjective, setObjectiveGradient, setHessian, petsc.TaoSetGradient
413
414        """
415        cdef PetscVec gvec = NULL
416        if g is not None: gvec = g.vec
417        if args is None: args = ()
418        if kargs is None: kargs = {}
419        context = (gradient, args, kargs)
420        self.set_attr("__gradient__", context)
421        CHKERR(TaoSetGradient(self.tao, gvec, TAO_Gradient, <void*>context))
422
423    def getObjective(self) -> TAOObjectiveFunction:
424        """Return the objective evaluation callback.
425
426        Not collective.
427
428        See Also
429        --------
430        setObjective, petsc.TaoGetObjective
431
432        """
433        cdef object objective = self.get_attr("__objective__")
434        return objective
435
436    def getGradient(self) -> tuple[Vec, TAOGradientFunction]:
437        """Return the vector used to store the gradient and the evaluation callback.
438
439        Not collective.
440
441        See Also
442        --------
443        setGradient, setHessian, petsc.TaoGetGradient
444
445        """
446        cdef Vec vec = Vec()
447        CHKERR(TaoGetGradient(self.tao, &vec.vec, NULL, NULL))
448        CHKERR(PetscINCREF(vec.obj))
449        cdef object gradient = self.get_attr("__gradient__")
450        return (vec, gradient)
451
452    def setObjectiveGradient(self, objgrad: TAOObjectiveGradientFunction, Vec g=None, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
453        """Set the objective function and gradient evaluation callback.
454
455        Logically collective.
456
457        Parameters
458        ----------
459        objgrad
460            The objective function and gradient callback.
461        g
462            The vector to store the gradient.
463        args
464            Positional arguments for the callback.
465        kargs
466            Keyword arguments for the callback.
467
468        See Also
469        --------
470        setObjective, setGradient, setHessian, getObjectiveAndGradient
471        petsc.TaoSetObjectiveAndGradient
472
473        """
474        cdef PetscVec gvec = NULL
475        if g is not None: gvec = g.vec
476        if args is None: args = ()
477        if kargs is None: kargs = {}
478        context = (objgrad, args, kargs)
479        self.set_attr("__objgrad__", context)
480        CHKERR(TaoSetObjectiveAndGradient(self.tao, gvec, TAO_ObjGrad, <void*>context))
481
482    def getObjectiveAndGradient(self) -> tuple[Vec, TAOObjectiveGradientFunction]:
483        """Return the vector used to store the gradient and the evaluation callback.
484
485        Not collective.
486
487        See Also
488        --------
489        setObjectiveGradient, petsc.TaoGetObjectiveAndGradient
490
491        """
492        cdef Vec vec = Vec()
493        CHKERR(TaoGetObjectiveAndGradient(self.tao, &vec.vec, NULL, NULL))
494        CHKERR(PetscINCREF(vec.obj))
495        cdef object objgrad = self.get_attr("__objgrad__")
496        return (vec, objgrad)
497
498    def setVariableBounds(self, varbounds: tuple[Vec, Vec] | TAOVariableBoundsFunction, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
499        """Set the lower and upper bounds for the optimization problem.
500
501        Logically collective.
502
503        Parameters
504        ----------
505        varbounds
506            Either a tuple of `Vec` or a `TAOVariableBoundsFunction` callback.
507        args
508            Positional arguments for the callback.
509        kargs
510            Keyword arguments for the callback.
511
512        See Also
513        --------
514        petsc.TaoSetVariableBounds, petsc.TaoSetVariableBoundsRoutine
515
516        """
517        cdef Vec xl = None, xu = None
518        if (isinstance(varbounds, list) or isinstance(varbounds, tuple)):
519            ol, ou = varbounds
520            xl = <Vec?> ol; xu = <Vec?> ou
521            CHKERR(TaoSetVariableBounds(self.tao, xl.vec, xu.vec))
522            return
523        if isinstance(varbounds, Vec): # FIXME
524            ol = varbounds; ou = args
525            xl = <Vec?> ol; xu = <Vec?> ou
526            CHKERR(TaoSetVariableBounds(self.tao, xl.vec, xu.vec))
527            return
528        if args is None: args = ()
529        if kargs is None: kargs = {}
530        context = (varbounds, args, kargs)
531        self.set_attr("__varbounds__", context)
532        CHKERR(TaoSetVariableBoundsRoutine(self.tao, TAO_VarBounds, <void*>context))
533
534    def setConstraints(self, constraints: TAOConstraintsFunction, Vec C=None, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
535        """Set the callback to compute constraints.
536
537        Logically collective.
538
539        Parameters
540        ----------
541        constraints
542            The callback.
543        C
544            The vector to hold the constraints.
545        args
546            Positional arguments for the callback.
547        kargs
548            Keyword arguments for the callback.
549
550        See Also
551        --------
552        petsc.TaoSetConstraintsRoutine
553
554        """
555        cdef PetscVec Cvec = NULL
556        if C is not None: Cvec = C.vec
557        if args is None: args = ()
558        if kargs is None: kargs = {}
559        context = (constraints, args, kargs)
560        self.set_attr("__constraints__", context)
561        CHKERR(TaoSetConstraintsRoutine(self.tao, Cvec, TAO_Constraints, <void*>context))
562
563    def setHessian(self, hessian: TAOHessianFunction, Mat H=None, Mat P=None,
564                   args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
565        """Set the callback to compute the Hessian matrix.
566
567        Logically collective.
568
569        Parameters
570        ----------
571        hessian
572            The Hessian callback.
573        H
574            The matrix to store the Hessian.
575        P
576            The matrix to construct the preconditioner.
577        args
578            Positional arguments for the callback.
579        kargs
580            Keyword arguments for the callback.
581
582        See Also
583        --------
584        getHessian, setObjective, setObjectiveGradient, setGradient
585        petsc.TaoSetHessian
586
587        """
588        cdef PetscMat Hmat = NULL
589        if H is not None: Hmat = H.mat
590        cdef PetscMat Pmat = Hmat
591        if P is not None: Pmat = P.mat
592        if args is None: args = ()
593        if kargs is None: kargs = {}
594        context = (hessian, args, kargs)
595        self.set_attr("__hessian__", context)
596        CHKERR(TaoSetHessian(self.tao, Hmat, Pmat, TAO_Hessian, <void*>context))
597
598    def getHessian(self) -> tuple[Mat, Mat, TAOHessianFunction]:
599        """Return the matrices used to store the Hessian and the evaluation callback.
600
601        Not collective.
602
603        See Also
604        --------
605        setHessian, petsc.TaoGetHessian
606
607        """
608        cdef Mat J = Mat()
609        cdef Mat P = Mat()
610        CHKERR(TaoGetHessian(self.tao, &J.mat, &P.mat, NULL, NULL))
611        CHKERR(PetscINCREF(J.obj))
612        CHKERR(PetscINCREF(P.obj))
613        cdef object hessian = self.get_attr("__hessian__")
614        return (J, P, hessian)
615
616    def setJacobian(self, jacobian: TAOJacobianFunction, Mat J=None, Mat P=None,
617                    args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
618        """Set the callback to compute the Jacobian.
619
620        Logically collective.
621
622        Parameters
623        ----------
624        jacobian
625            The Jacobian callback.
626        J
627            The matrix to store the Jacobian.
628        P
629            The matrix to construct the preconditioner.
630        args
631            Positional arguments for the callback.
632        kargs
633            Keyword arguments for the callback.
634
635        See Also
636        --------
637        petsc.TaoSetJacobianRoutine
638
639        """
640        cdef PetscMat Jmat = NULL
641        if J is not None: Jmat = J.mat
642        cdef PetscMat Pmat = Jmat
643        if P is not None: Pmat = P.mat
644        if args is None: args = ()
645        if kargs is None: kargs = {}
646        context = (jacobian, args, kargs)
647        self.set_attr("__jacobian__", context)
648        CHKERR(TaoSetJacobianRoutine(self.tao, Jmat, Pmat, TAO_Jacobian, <void*>context))
649
650    def setStateDesignIS(self, IS state=None, IS design=None) -> None:
651        """Set the index sets indicating state and design variables.
652
653        Collective.
654
655        See Also
656        --------
657        petsc.TaoSetStateDesignIS
658
659        """
660        cdef PetscIS s_is = NULL, d_is = NULL
661        if state  is not None: s_is = state.iset
662        if design is not None: d_is = design.iset
663        CHKERR(TaoSetStateDesignIS(self.tao, s_is, d_is))
664
665    def setJacobianState(self, jacobian_state, Mat J=None, Mat P=None, Mat I=None,
666                         args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
667        """Set Jacobian state callback.
668
669        Logically collective.
670
671        See Also
672        --------
673        petsc.TaoSetJacobianStateRoutine
674
675        """
676        cdef PetscMat Jmat = NULL
677        if J is not None: Jmat = J.mat
678        cdef PetscMat Pmat = Jmat
679        if P is not None: Pmat = P.mat
680        cdef PetscMat Imat = NULL
681        if I is not None: Imat = I.mat
682        if args is None: args = ()
683        if kargs is None: kargs = {}
684        context = (jacobian_state, args, kargs)
685        self.set_attr("__jacobian_state__", context)
686        CHKERR(TaoSetJacobianStateRoutine(self.tao, Jmat, Pmat, Imat,
687                                          TAO_JacobianState, <void*>context))
688
689    def setJacobianDesign(self, jacobian_design, Mat J=None,
690                          args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
691        """Set Jacobian design callback.
692
693        Logically collective.
694
695        See Also
696        --------
697        petsc.TaoSetJacobianDesignRoutine
698
699        """
700        cdef PetscMat Jmat = NULL
701        if J is not None: Jmat = J.mat
702        if args is None: args = ()
703        if kargs is None: kargs = {}
704        context = (jacobian_design, args, kargs)
705        self.set_attr("__jacobian_design__", context)
706        CHKERR(TaoSetJacobianDesignRoutine(self.tao, Jmat,
707                                           TAO_JacobianDesign, <void*>context))
708
709    def getLMVMMat(self) -> Mat:
710        """Get the LMVM matrix.
711
712        Not collective.
713
714        See Also
715        --------
716        setLMVMMat, petsc.TaoGetLMVMMatrix
717
718        """
719        cdef Mat M = Mat()
720        CHKERR(TaoGetLMVMMatrix(self.tao, &M.mat))
721        CHKERR(PetscINCREF(M.obj))
722        return M
723
724    def setLMVMMat(self, Mat M) -> None:
725        """Set the LMVM matrix.
726
727        Logically collective.
728
729        See Also
730        --------
731        getLMVMMat, petsc.TaoSetLMVMMatrix
732
733        """
734        cdef PetscMat ctype = M.mat
735        CHKERR(TaoSetLMVMMatrix(self.tao, ctype))
736
737    def setEqualityConstraints(self, equality_constraints, Vec c,
738                               args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
739        """Set equality constraints callback.
740
741        Logically collective.
742
743        See Also
744        --------
745        petsc.TaoSetEqualityConstraintsRoutine
746
747        """
748        if args is None: args = ()
749        if kargs is None: kargs = {}
750        context = (equality_constraints, args, kargs)
751        self.set_attr("__equality_constraints__", context)
752        CHKERR(TaoSetEqualityConstraintsRoutine(self.tao, c.vec,
753                                                TAO_EqualityConstraints, <void*>context))
754
755    def getEqualityConstraints(self) -> tuple[Vec, tuple[TAOConstraintsFunction, tuple[Any, ...] | None, dict[str, Any] | None]]:
756        """Return tuple holding vector and callback of equality constraints.
757
758        Not collective.
759
760        See Also
761        --------
762        setEqualityConstraints, petsc.TaoGetEqualityConstraintsRoutine
763        """
764        cdef Vec c = Vec()
765        CHKERR(TaoGetEqualityConstraintsRoutine(self.tao, &c.vec, NULL, NULL))
766        CHKERR(PetscINCREF(c.obj))
767        cdef object equality_constraints = self.get_attr("__equality_constraints__")
768        return (c, equality_constraints)
769
770    def setJacobianEquality(self, jacobian_equality, Mat J=None, Mat P=None,
771                            args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
772        """Set Jacobian equality constraints callback.
773
774        Logically collective.
775
776        See Also
777        --------
778        petsc.TaoSetJacobianEqualityRoutine
779
780        """
781        cdef PetscMat Jmat = NULL
782        if J is not None: Jmat = J.mat
783        cdef PetscMat Pmat = Jmat
784        if P is not None: Pmat = P.mat
785        if args is None: args = ()
786        if kargs is None: kargs = {}
787        context = (jacobian_equality, args, kargs)
788        self.set_attr("__jacobian_equality__", context)
789        CHKERR(TaoSetJacobianEqualityRoutine(self.tao, Jmat, Pmat,
790                                             TAO_JacobianEquality, <void*>context))
791
792    def getJacobianEquality(self) -> tuple[Mat, Mat, tuple[TAOConstraintsJacobianFunction, tuple[Any, ...] | None,
793                                           dict[str, Any] | None]]:
794        """Return matrix, precon matrix and callback of equality constraints Jacobian.
795
796        Not collective.
797
798        See Also
799        --------
800        setJacobianEquality, petsc.TaoGetJacobianEqualityRoutine
801        """
802        cdef Mat J = Mat()
803        cdef Mat Jpre = Mat()
804        CHKERR(TaoGetJacobianEqualityRoutine(self.tao, &J.mat, &Jpre.mat, NULL, NULL))
805        CHKERR(PetscINCREF(J.obj))
806        CHKERR(PetscINCREF(Jpre.obj))
807        cdef object jacobian_equality = self.get_attr("__jacobian_equality__")
808        return (J, Jpre, jacobian_equality)
809
810    def setInequalityConstraints(self, inequality_constraints, Vec c,
811                                 args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
812        """Set inequality constraints callback.
813
814        Logically collective.
815
816        See Also
817        --------
818        petsc.TaoSetInequalityConstraintsRoutine
819
820        """
821        if args is None: args = ()
822        if kargs is None: kargs = {}
823        context = (inequality_constraints, args, kargs)
824        self.set_attr("__inequality_constraints__", context)
825        CHKERR(TaoSetInequalityConstraintsRoutine(self.tao, c.vec,
826                                                  TAO_InequalityConstraints, <void*>context))
827
828    def getInequalityConstraints(self) -> tuple[Vec, tuple[TAOConstraintsFunction, tuple[Any, ...] | None, dict[str, Any] | None]]:
829        """Return tuple holding vector and callback of inequality constraints.
830
831        Not collective.
832
833        See Also
834        --------
835        setInequalityConstraints, petsc.TaoGetInequalityConstraintsRoutine
836        """
837        cdef Vec c = Vec()
838        CHKERR(TaoGetInequalityConstraintsRoutine(self.tao, &c.vec, NULL, NULL))
839        CHKERR(PetscINCREF(c.obj))
840        cdef object inequality_constraints = self.get_attr("__inequality_constraints__")
841        return (c, inequality_constraints)
842
843    def setJacobianInequality(self, jacobian_inequality, Mat J=None, Mat P=None,
844                              args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
845        """Set Jacobian inequality constraints callback.
846
847        Logically collective.
848
849        See Also
850        --------
851        petsc.TaoSetJacobianInequalityRoutine
852
853        """
854        cdef PetscMat Jmat = NULL
855        if J is not None: Jmat = J.mat
856        cdef PetscMat Pmat = Jmat
857        if P is not None: Pmat = P.mat
858        if args is None: args = ()
859        if kargs is None: kargs = {}
860        context = (jacobian_inequality, args, kargs)
861        self.set_attr("__jacobian_inequality__", context)
862        CHKERR(TaoSetJacobianInequalityRoutine(self.tao, Jmat, Pmat,
863                                               TAO_JacobianInequality, <void*>context))
864
865    def getJacobianInequality(self) -> tuple[Mat, Mat, tuple[TAOConstraintsJacobianFunction, tuple[Any, ...] | None, dict[str, Any] | None]]:
866        """Return matrix, precon matrix and callback of ineq. constraints Jacobian.
867
868        Not collective.
869
870        See Also
871        --------
872        setJacobianInequality, petsc.TaoGetJacobianInequalityRoutine
873        """
874        cdef Mat J = Mat()
875        cdef Mat Jpre = Mat()
876        CHKERR(TaoGetJacobianInequalityRoutine(self.tao, &J.mat, &Jpre.mat, NULL, NULL))
877        CHKERR(PetscINCREF(J.obj))
878        CHKERR(PetscINCREF(Jpre.obj))
879        cdef object jacobian_inequality = self.get_attr("__jacobian_inequality__")
880        return (J, Jpre, jacobian_inequality)
881
882    def setUpdate(self, update: TAOUpdateFunction, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
883        """Set the callback to compute update at each optimization step.
884
885        Logically collective.
886
887        Parameters
888        ----------
889        update
890            The update callback or `None` to reset it.
891        args
892            Positional arguments for the callback.
893        kargs
894            Keyword arguments for the callback.
895
896        See Also
897        --------
898        getUpdate, petsc.TaoSetUpdate
899
900        """
901        if update is not None:
902            if args  is None: args  = ()
903            if kargs is None: kargs = {}
904            context = (update, args, kargs)
905            self.set_attr('__update__', context)
906            CHKERR(TaoSetUpdate(self.tao, TAO_Update, <void*>context))
907        else:
908            self.set_attr('__update__', None)
909            CHKERR(TaoSetUpdate(self.tao, NULL, NULL))
910
911    def getUpdate(self) -> tuple[TAOUpdateFunction, tuple[Any, ...], dict[str, Any]]:
912        """Return the callback to compute the update.
913
914        Not collective.
915
916        See Also
917        --------
918        setUpdate
919
920        """
921        return self.get_attr('__update__')
922
923    # --------------
924
925    def computeObjective(self, Vec x) -> float:
926        """Compute the value of the objective function.
927
928        Collective.
929
930        Parameters
931        ----------
932        x
933            The parameter vector.
934
935        See Also
936        --------
937        setObjective, petsc.TaoComputeObjective
938
939        """
940        cdef PetscReal f = 0
941        CHKERR(TaoComputeObjective(self.tao, x.vec, &f))
942        return toReal(f)
943
944    def computeResidual(self, Vec x, Vec f) -> None:
945        """Compute the residual.
946
947        Collective.
948
949        Parameters
950        ----------
951        x
952            The parameter vector.
953        f
954            The output vector.
955
956        See Also
957        --------
958        setResidual, petsc.TaoComputeResidual
959
960        """
961        CHKERR(TaoComputeResidual(self.tao, x.vec, f.vec))
962
963    def computeGradient(self, Vec x, Vec g) -> None:
964        """Compute the gradient of the objective function.
965
966        Collective.
967
968        Parameters
969        ----------
970        x
971            The parameter vector.
972        g
973            The output gradient vector.
974
975        See Also
976        --------
977        setGradient, petsc.TaoComputeGradient
978
979        """
980        CHKERR(TaoComputeGradient(self.tao, x.vec, g.vec))
981
982    def computeObjectiveGradient(self, Vec x, Vec g) -> float:
983        """Compute the gradient of the objective function and its value.
984
985        Collective.
986
987        Parameters
988        ----------
989        x
990            The parameter vector.
991        g
992            The output gradient vector.
993
994        See Also
995        --------
996        setObjectiveGradient, setGradient, setObjective
997        petsc.TaoComputeObjectiveAndGradient
998
999        """
1000        cdef PetscReal f = 0
1001        CHKERR(TaoComputeObjectiveAndGradient(self.tao, x.vec, &f, g.vec))
1002        return toReal(f)
1003
1004    def computeDualVariables(self, Vec xl, Vec xu) -> None:
1005        """Compute the dual vectors corresponding to variables' bounds.
1006
1007        Collective.
1008
1009        See Also
1010        --------
1011        petsc.TaoComputeDualVariables
1012
1013        """
1014        CHKERR(TaoComputeDualVariables(self.tao, xl.vec, xu.vec))
1015
1016    def computeVariableBounds(self, Vec xl, Vec xu) -> None:
1017        """Compute the vectors corresponding to variables' bounds.
1018
1019        Collective.
1020
1021        See Also
1022        --------
1023        setVariableBounds, petsc.TaoComputeVariableBounds
1024
1025        """
1026        CHKERR(TaoComputeVariableBounds(self.tao))
1027        cdef PetscVec Lvec = NULL, Uvec = NULL
1028        CHKERR(TaoGetVariableBounds(self.tao, &Lvec, &Uvec))
1029        if xl.vec != NULL:
1030            if Lvec != NULL:
1031                CHKERR(VecCopy(Lvec, xl.vec))
1032            else:
1033                CHKERR(VecSet(xl.vec, <PetscScalar>PETSC_NINFINITY))
1034        if xu.vec != NULL:
1035            if Uvec != NULL:
1036                CHKERR(VecCopy(Uvec, xu.vec))
1037            else:
1038                CHKERR(VecSet(xu.vec, <PetscScalar>PETSC_INFINITY))
1039
1040    def computeConstraints(self, Vec x, Vec c) -> None:
1041        """Compute the vector corresponding to the constraints.
1042
1043        Collective.
1044
1045        Parameters
1046        ----------
1047        x
1048            The parameter vector.
1049        c
1050            The output constraints vector.
1051
1052        See Also
1053        --------
1054        setVariableBounds, petsc.TaoComputeVariableBounds
1055
1056        """
1057        CHKERR(TaoComputeConstraints(self.tao, x.vec, c.vec))
1058
1059    def computeHessian(self, Vec x, Mat H, Mat P=None) -> None:
1060        """Compute the Hessian of the objective function.
1061
1062        Collective.
1063
1064        Parameters
1065        ----------
1066        x
1067            The parameter vector.
1068        H
1069            The output Hessian matrix.
1070        P
1071            The output Hessian matrix used to construct the preconditioner.
1072
1073        See Also
1074        --------
1075        setHessian, petsc.TaoComputeHessian
1076
1077        """
1078        cdef PetscMat hmat = H.mat, pmat = H.mat
1079        if P is not None: pmat = P.mat
1080        CHKERR(TaoComputeHessian(self.tao, x.vec, hmat, pmat))
1081
1082    def computeJacobian(self, Vec x, Mat J, Mat P=None) -> None:
1083        """Compute the Jacobian.
1084
1085        Collective.
1086
1087        Parameters
1088        ----------
1089        x
1090            The parameter vector.
1091        J
1092            The output Jacobian matrix.
1093        P
1094            The output Jacobian matrix used to construct the preconditioner.
1095
1096        See Also
1097        --------
1098        setJacobian, petsc.TaoComputeJacobian
1099
1100        """
1101        cdef PetscMat jmat = J.mat, pmat = J.mat
1102        if P is not None: pmat = P.mat
1103        CHKERR(TaoComputeJacobian(self.tao, x.vec, jmat, pmat))
1104
1105    # --------------
1106
1107    def setTolerances(self, gatol: float | None = None, grtol: float | None = None, gttol: float | None = None) -> None:
1108        """Set the tolerance parameters used in the solver convergence tests.
1109
1110        Collective.
1111
1112        Parameters
1113        ----------
1114        gatol
1115            The absolute norm of the gradient, or `DETERMINE`
1116            to use the value when the object's type was set.
1117            Defaults to `CURRENT`.
1118        grtol
1119            The relative norm of the gradient with respect to the
1120            initial norm of the objective, or `DETERMINE` to
1121            use the value when the object's type was set. Defaults
1122            to `CURRENT`.
1123        gttol
1124            The relative norm of the gradient with respect to the
1125            initial norm of the gradient, or `DETERMINE` to
1126            use the value when the object's type was set. Defaults
1127            to `CURRENT`.
1128
1129        See Also
1130        --------
1131        getTolerances, petsc.TaoSetTolerances
1132
1133        """
1134        cdef PetscReal _gatol=PETSC_CURRENT, _grtol=PETSC_CURRENT, _gttol=PETSC_CURRENT
1135        if gatol is not None: _gatol = asReal(gatol)
1136        if grtol is not None: _grtol = asReal(grtol)
1137        if gttol is not None: _gttol = asReal(gttol)
1138        CHKERR(TaoSetTolerances(self.tao, _gatol, _grtol, _gttol))
1139
1140    def getTolerances(self) -> tuple[float, float, float]:
1141        """Return the tolerance parameters used in the solver convergence tests.
1142
1143        Not collective.
1144
1145        Returns
1146        -------
1147        gatol : float
1148            The absolute norm of the gradient.
1149        grtol : float
1150            The relative norm of the gradient.
1151        gttol : float
1152            The relative norm of the gradient with respect to the
1153            initial norm of the gradient.
1154
1155        See Also
1156        --------
1157        setTolerances, petsc.TaoGetTolerances
1158
1159        """
1160        cdef PetscReal _gatol=PETSC_CURRENT, _grtol=PETSC_CURRENT, _gttol=PETSC_CURRENT
1161        CHKERR(TaoGetTolerances(self.tao, &_gatol, &_grtol, &_gttol))
1162        return (toReal(_gatol), toReal(_grtol), toReal(_gttol))
1163
1164    def setMaximumIterations(self, mit: int) -> float:
1165        """Set the maximum number of solver iterations.
1166
1167        Collective.
1168
1169        See Also
1170        --------
1171        setTolerances, petsc.TaoSetMaximumIterations
1172
1173        """
1174        cdef PetscInt _mit = asInt(mit)
1175        CHKERR(TaoSetMaximumIterations(self.tao, _mit))
1176
1177    def getMaximumIterations(self) -> int:
1178        """Return the maximum number of solver iterations.
1179
1180        Not collective.
1181
1182        See Also
1183        --------
1184        setMaximumIterations, petsc.TaoGetMaximumIterations
1185
1186        """
1187        cdef PetscInt _mit = PETSC_DEFAULT
1188        CHKERR(TaoGetMaximumIterations(self.tao, &_mit))
1189        return toInt(_mit)
1190
1191    def setMaximumFunctionEvaluations(self, mit: int) -> None:
1192        """Set the maximum number of objective evaluations within the solver.
1193
1194        Collective.
1195
1196        See Also
1197        --------
1198        setMaximumIterations, petsc.TaoSetMaximumFunctionEvaluations
1199
1200        """
1201        cdef PetscInt _mit = asInt(mit)
1202        CHKERR(TaoSetMaximumFunctionEvaluations(self.tao, _mit))
1203
1204    def getMaximumFunctionEvaluations(self) -> int:
1205        """Return the maximum number of objective evaluations within the solver.
1206
1207        Not collective.
1208
1209        See Also
1210        --------
1211        setMaximumFunctionEvaluations, petsc.TaoGetMaximumFunctionEvaluations
1212
1213        """
1214        cdef PetscInt _mit = PETSC_DEFAULT
1215        CHKERR(TaoGetMaximumFunctionEvaluations(self.tao, &_mit))
1216        return toInt(_mit)
1217
1218    def setConstraintTolerances(self, catol: float | None = None, crtol: float | None = None) -> None:
1219        """Set the constraints tolerance parameters used in the solver convergence tests.
1220
1221        Collective.
1222
1223        Parameters
1224        ----------
1225        catol
1226            The absolute norm of the constraints, or `DETERMINE`
1227            to use the value when the object's type was set. Defaults
1228            to `CURRENT`.
1229        crtol
1230            The relative norm of the constraints, or `DETERMINE`
1231            to use the value when the object's type was set. Defaults
1232            to `CURRENT`.
1233
1234        See Also
1235        --------
1236        getConstraintTolerances, petsc.TaoSetConstraintTolerances
1237
1238        """
1239        cdef PetscReal _catol=PETSC_CURRENT, _crtol=PETSC_CURRENT
1240        if catol is not None: _catol = asReal(catol)
1241        if crtol is not None: _crtol = asReal(crtol)
1242        CHKERR(TaoSetConstraintTolerances(self.tao, _catol, _crtol))
1243
1244    def getConstraintTolerances(self) -> tuple[float, float]:
1245        """Return the constraints tolerance parameters used in the convergence tests.
1246
1247        Not collective.
1248
1249        Returns
1250        -------
1251        catol : float
1252            The absolute norm of the constraints.
1253        crtol : float
1254            The relative norm of the constraints.
1255
1256        See Also
1257        --------
1258        setConstraintTolerances, petsc.TaoGetConstraintTolerances
1259
1260        """
1261        cdef PetscReal _catol=PETSC_DEFAULT, _crtol=PETSC_DEFAULT
1262        CHKERR(TaoGetConstraintTolerances(self.tao, &_catol, &_crtol))
1263        return (toReal(_catol), toReal(_crtol))
1264
1265    def setConvergenceTest(self, converged: TAOConvergedFunction | None, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
1266        """Set the callback used to test for solver convergence.
1267
1268        Logically collective.
1269
1270        Parameters
1271        ----------
1272        converged
1273            The callback. If `None`, reset to the default convergence test.
1274        args
1275            Positional arguments for the callback.
1276        kargs
1277            Keyword arguments for the callback.
1278
1279        See Also
1280        --------
1281        getConvergenceTest, petsc.TaoSetConvergenceTest
1282
1283        """
1284        if converged is None:
1285            CHKERR(TaoSetConvergenceTest(self.tao, TaoDefaultConvergenceTest, NULL))
1286            self.set_attr('__converged__', None)
1287        else:
1288            if args is None: args = ()
1289            if kargs is None: kargs = {}
1290            self.set_attr('__converged__', (converged, args, kargs))
1291            CHKERR(TaoSetConvergenceTest(self.tao, TAO_Converged, NULL))
1292
1293    def getConvergenceTest(self) -> tuple[TAOConvergedFunction, tuple[Any, ...], dict[str, Any]]:
1294        """Return the callback used to test for solver convergence.
1295
1296        Not collective.
1297
1298        See Also
1299        --------
1300        setConvergenceTest
1301
1302        """
1303        return self.get_attr('__converged__')
1304
1305    def setConvergedReason(self, reason: ConvergedReason) -> None:
1306        """Set the termination flag.
1307
1308        Collective.
1309
1310        See Also
1311        --------
1312        getConvergedReason, petsc.TaoSetConvergedReason
1313
1314        """
1315        cdef PetscTAOConvergedReason creason = reason
1316        CHKERR(TaoSetConvergedReason(self.tao, creason))
1317
1318    def getConvergedReason(self) -> ConvergedReason:
1319        """Return the termination flag.
1320
1321        Not collective.
1322
1323        See Also
1324        --------
1325        setConvergedReason, petsc.TaoGetConvergedReason
1326
1327        """
1328        cdef PetscTAOConvergedReason creason = TAO_CONTINUE_ITERATING
1329        CHKERR(TaoGetConvergedReason(self.tao, &creason))
1330        return creason
1331
1332    def setMonitor(self, monitor: TAOMonitorFunction, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
1333        """Set the callback used to monitor solver convergence.
1334
1335        Logically collective.
1336
1337        Parameters
1338        ----------
1339        monitor
1340            The callback.
1341        args
1342            Positional arguments for the callback.
1343        kargs
1344            Keyword arguments for the callback.
1345
1346        See Also
1347        --------
1348        getMonitor, petsc.TaoMonitorSet
1349
1350        """
1351        if monitor is None: return
1352        cdef object monitorlist = self.get_attr('__monitor__')
1353        if args  is None: args  = ()
1354        if kargs is None: kargs = {}
1355        if monitorlist is None:
1356            CHKERR(TaoMonitorSet(self.tao, TAO_Monitor, NULL, NULL))
1357            self.set_attr('__monitor__',  [(monitor, args, kargs)])
1358        else:
1359            monitorlist.append((monitor, args, kargs))
1360
1361    def getMonitor(self) -> list[tuple[TAOMonitorFunction, tuple[Any, ...], dict[str, Any]]]:
1362        """Return the callback used to monitor solver convergence.
1363
1364        Not collective.
1365
1366        See Also
1367        --------
1368        setMonitor
1369
1370        """
1371        return self.get_attr('__monitor__')
1372
1373    def cancelMonitor(self) -> None:
1374        """Cancel all the monitors of the solver.
1375
1376        Logically collective.
1377
1378        See Also
1379        --------
1380        setMonitor, petsc.TaoMonitorCancel
1381
1382        """
1383        CHKERR(TaoMonitorCancel(self.tao))
1384        self.set_attr('__monitor__',  None)
1385
1386    # Tao overwrites these statistics. Copy user defined only if present
1387    def monitor(self, its: int | None = None, f: float | None = None, res: float | None = None, cnorm: float | None = None, step: float | None = None) -> None:
1388        """Monitor the solver.
1389
1390        Collective.
1391
1392        This function should be called without arguments,
1393        unless users want to modify the values internally stored by the solver.
1394
1395        Parameters
1396        ----------
1397        its
1398            Current number of iterations
1399            or `None` to use the value stored internally by the solver.
1400        f
1401            Current value of the objective function
1402            or `None` to use the value stored internally by the solver.
1403        res
1404            Current value of the residual norm
1405            or `None` to use the value stored internally by the solver.
1406        cnorm
1407            Current value of the constrains norm
1408            or `None` to use the value stored internally by the solver.
1409        step
1410            Current value of the step
1411            or `None` to use the value stored internally by the solver.
1412
1413        See Also
1414        --------
1415        setMonitor, petsc.TaoMonitor
1416
1417        """
1418        cdef PetscInt cits = 0
1419        cdef PetscReal cf = 0.0
1420        cdef PetscReal cres = 0.0
1421        cdef PetscReal ccnorm = 0.0
1422        cdef PetscReal cstep = 0.0
1423        CHKERR(TaoGetSolutionStatus(self.tao, &cits, &cf, &cres, &ccnorm, &cstep, NULL))
1424        if its is not None:
1425            cits = asInt(its)
1426        if f is not None:
1427            cf = asReal(f)
1428        if res is not None:
1429            cres = asReal(res)
1430        if cnorm is not None:
1431            ccnorm = asReal(cnorm)
1432        if step is not None:
1433            cstep = asReal(step)
1434        CHKERR(TaoMonitor(self.tao, cits, cf, cres, ccnorm, cstep))
1435
1436    #
1437
1438    def solve(self, Vec x=None) -> None:
1439        """Solve the optimization problem.
1440
1441        Collective.
1442
1443        Parameters
1444        ----------
1445        x
1446            The starting vector or `None` to use the vector stored internally.
1447
1448        See Also
1449        --------
1450        setSolution, getSolution, petsc.TaoSolve
1451
1452        """
1453        if x is not None:
1454            CHKERR(TaoSetSolution(self.tao, x.vec))
1455        CHKERR(TaoSolve(self.tao))
1456
1457    def getSolution(self) -> Vec:
1458        """Return the vector holding the solution.
1459
1460        Not collective.
1461
1462        See Also
1463        --------
1464        setSolution, petsc.TaoGetSolution
1465
1466        """
1467        cdef Vec vec = Vec()
1468        CHKERR(TaoGetSolution(self.tao, &vec.vec))
1469        CHKERR(PetscINCREF(vec.obj))
1470        return vec
1471
1472    def setGradientNorm(self, Mat mat) -> None:
1473        """Set the matrix used to compute inner products.
1474
1475        Collective.
1476
1477        See Also
1478        --------
1479        getGradientNorm, petsc.TaoSetGradientNorm
1480
1481        """
1482        CHKERR(TaoSetGradientNorm(self.tao, mat.mat))
1483
1484    def getGradientNorm(self) -> Mat:
1485        """Return the matrix used to compute inner products.
1486
1487        Not collective.
1488
1489        See Also
1490        --------
1491        setGradientNorm, petsc.TaoGetGradientNorm
1492
1493        """
1494        cdef Mat mat = Mat()
1495        CHKERR(TaoGetGradientNorm(self.tao, &mat.mat))
1496        CHKERR(PetscINCREF(mat.obj))
1497        return mat
1498
1499    def setLMVMH0(self, Mat mat) -> None:
1500        """Set the initial Hessian for the quasi-Newton approximation.
1501
1502        Collective.
1503
1504        See Also
1505        --------
1506        getLMVMH0, petsc.TaoLMVMSetH0
1507
1508        """
1509        CHKERR(TaoLMVMSetH0(self.tao, mat.mat))
1510
1511    def getLMVMH0(self) -> Mat:
1512        """Return the initial Hessian for the quasi-Newton approximation.
1513
1514        Not collective.
1515
1516        See Also
1517        --------
1518        setLMVMH0, petsc.TaoLMVMGetH0
1519
1520        """
1521        cdef Mat mat = Mat()
1522        CHKERR(TaoLMVMGetH0(self.tao, &mat.mat))
1523        CHKERR(PetscINCREF(mat.obj))
1524        return mat
1525
1526    def getLMVMH0KSP(self) -> KSP:
1527        """Return the `KSP` for the inverse of the initial Hessian approximation.
1528
1529        Not collective.
1530
1531        See Also
1532        --------
1533        setLMVMH0, petsc.TaoLMVMGetH0KSP
1534
1535        """
1536        cdef KSP ksp = KSP()
1537        CHKERR(TaoLMVMGetH0KSP(self.tao, &ksp.ksp))
1538        CHKERR(PetscINCREF(ksp.obj))
1539        return ksp
1540
1541    def getVariableBounds(self) -> tuple[Vec, Vec]:
1542        """Return the lower and upper bounds vectors.
1543
1544        Not collective.
1545
1546        See Also
1547        --------
1548        setVariableBounds, petsc.TaoGetVariableBounds
1549
1550        """
1551        cdef Vec xl = Vec(), xu = Vec()
1552        CHKERR(TaoGetVariableBounds(self.tao, &xl.vec, &xu.vec))
1553        CHKERR(PetscINCREF(xl.obj)); CHKERR(PetscINCREF(xu.obj))
1554        return (xl, xu)
1555
1556    def setBNCGType(self, cg_type: BNCGType) -> None:
1557        """Set the type of the BNCG solver.
1558
1559        Collective.
1560
1561        See Also
1562        --------
1563        getBNCGType, petsc.TaoBNCGSetType
1564
1565        """
1566        cdef PetscTAOBNCGType ctype = cg_type
1567        CHKERR(TaoBNCGSetType(self.tao, ctype))
1568
1569    def getBNCGType(self) -> BNCGType:
1570        """Return the type of the BNCG solver.
1571
1572        Not collective.
1573
1574        See Also
1575        --------
1576        setBNCGType, petsc.TaoBNCGGetType
1577
1578        """
1579        cdef PetscTAOBNCGType cg_type = TAO_BNCG_SSML_BFGS
1580        CHKERR(TaoBNCGGetType(self.tao, &cg_type))
1581        return cg_type
1582
1583    def setIterationNumber(self, its: int) -> None:
1584        """Set the current iteration number.
1585
1586        Collective.
1587
1588        See Also
1589        --------
1590        getIterationNumber, petsc.TaoSetIterationNumber
1591
1592        """
1593        cdef PetscInt ival = asInt(its)
1594        CHKERR(TaoSetIterationNumber(self.tao, ival))
1595
1596    def getIterationNumber(self) -> int:
1597        """Return the current iteration number.
1598
1599        Not collective.
1600
1601        See Also
1602        --------
1603        setIterationNumber, petsc.TaoGetIterationNumber
1604
1605        """
1606        cdef PetscInt its=0
1607        CHKERR(TaoGetIterationNumber(self.tao, &its))
1608        return toInt(its)
1609
1610    def getObjectiveValue(self) -> float:
1611        """Return the current value of the objective function.
1612
1613        Not collective.
1614
1615        See Also
1616        --------
1617        setObjective, petsc.TaoGetSolutionStatus
1618
1619        """
1620        cdef PetscReal fval=0
1621        CHKERR(TaoGetSolutionStatus(self.tao, NULL, &fval, NULL, NULL, NULL, NULL))
1622        return toReal(fval)
1623
1624    getFunctionValue = getObjectiveValue
1625
1626    def getConvergedReason(self) -> ConvergedReason:
1627        """Return the reason for the solver convergence.
1628
1629        Not collective.
1630
1631        See Also
1632        --------
1633        petsc.TaoGetConvergedReason
1634
1635        """
1636        cdef PetscTAOConvergedReason reason = TAO_CONTINUE_ITERATING
1637        CHKERR(TaoGetConvergedReason(self.tao, &reason))
1638        return reason
1639
1640    def getSolutionNorm(self) -> tuple[float, float, float]:
1641        """Return the objective function value and the norms of gradient and constraints.
1642
1643        Not collective.
1644
1645        Returns
1646        -------
1647        f : float
1648            Current value of the objective function.
1649        res : float
1650            Current value of the residual norm.
1651        cnorm : float
1652            Current value of the constrains norm.
1653
1654        See Also
1655        --------
1656        getSolutionStatus, petsc.TaoGetSolutionStatus
1657
1658        """
1659        cdef PetscReal gnorm=0
1660        cdef PetscReal cnorm=0
1661        cdef PetscReal fval=0
1662        CHKERR(TaoGetSolutionStatus(self.tao, NULL, &fval, &gnorm, &cnorm, NULL, NULL))
1663        return (toReal(fval), toReal(gnorm), toReal(cnorm))
1664
1665    def getSolutionStatus(self) -> tuple[int, float, float, float, float, ConvergedReason]:
1666        """Return the solution status.
1667
1668        Not collective.
1669
1670        Returns
1671        -------
1672        its : int
1673            Current number of iterations.
1674        f : float
1675            Current value of the objective function.
1676        res : float
1677            Current value of the residual norm.
1678        cnorm : float
1679            Current value of the constrains norm.
1680        step : float
1681            Current value of the step.
1682        reason : ConvergedReason
1683            Current value of converged reason.
1684
1685        See Also
1686        --------
1687        petsc.TaoGetSolutionStatus
1688
1689        """
1690        cdef PetscInt its=0
1691        cdef PetscReal fval=0, gnorm=0, cnorm=0, xdiff=0
1692        cdef PetscTAOConvergedReason reason = TAO_CONTINUE_ITERATING
1693        CHKERR(TaoGetSolutionStatus(self.tao, &its,
1694                                    &fval, &gnorm, &cnorm, &xdiff,
1695                                    &reason))
1696        return (toInt(its), toReal(fval),
1697                toReal(gnorm), toReal(cnorm),
1698                toReal(xdiff), reason)
1699
1700    def checkConverged(self) -> ConvergedReason:
1701        """Run convergence test and return converged reason.
1702
1703        Collective.
1704
1705        See Also
1706        --------
1707        converged
1708
1709        """
1710        cdef PetscTAOConvergedReason reason = TAO_CONTINUE_ITERATING
1711        CHKERR(TaoConverged(self.tao, &reason))
1712        return reason
1713
1714    def getKSP(self) -> KSP:
1715        """Return the linear solver used by the nonlinear solver.
1716
1717        Not collective.
1718
1719        See Also
1720        --------
1721        petsc.TaoGetKSP
1722
1723        """
1724        cdef KSP ksp = KSP()
1725        CHKERR(TaoGetKSP(self.tao, &ksp.ksp))
1726        CHKERR(PetscINCREF(ksp.obj))
1727        return ksp
1728
1729    # ALMM routines
1730
1731    def getALMMSubsolver(self) -> TAO:
1732        """Return the subsolver inside the ALMM solver.
1733
1734        Not collective.
1735
1736        See Also
1737        --------
1738        setALMMSubsolver, petsc.TaoALMMGetSubsolver
1739
1740        """
1741        cdef TAO subsolver = TAO()
1742        CHKERR(TaoALMMGetSubsolver(self.tao, &subsolver.tao))
1743        CHKERR(PetscINCREF(subsolver.obj))
1744        return subsolver
1745
1746    def getALMMType(self) -> ALMMType:
1747        """Return the type of the ALMM solver.
1748
1749        Not collective.
1750
1751        See Also
1752        --------
1753        setALMMType, petsc.TaoALMMGetType
1754
1755        """
1756        cdef PetscTAOALMMType almm_type = TAO_ALMM_PHR
1757        CHKERR(TaoALMMGetType(self.tao, &almm_type))
1758        return almm_type
1759
1760    def setALMMSubsolver(self, subsolver: TAO) -> None:
1761        """Set the subsolver inside the ALMM solver.
1762
1763        Logically collective.
1764
1765        See Also
1766        --------
1767        getALMMSubsolver, petsc.TaoALMMSetSubsolver
1768
1769        """
1770        cdef TAO ctype = subsolver
1771        CHKERR(TaoALMMSetSubsolver(self.tao, ctype.tao))
1772
1773    def setALMMType(self, tao_almm_type: ALMMType) -> None:
1774        """Set the ALMM type of the solver.
1775
1776        Logically collective.
1777
1778        Parameters
1779        ----------
1780        tao_almm_type
1781            The type of the solver.
1782
1783        See Also
1784        --------
1785        getALMMType, petsc.TaoALMMSetType
1786
1787        """
1788        cdef PetscTAOALMMType ctype = tao_almm_type
1789        CHKERR(TaoALMMSetType(self.tao, ctype))
1790
1791    # BRGN routines
1792
1793    def getBRGNSubsolver(self) -> TAO:
1794        """Return the subsolver inside the BRGN solver.
1795
1796        Not collective.
1797
1798        See Also
1799        --------
1800        petsc.TaoBRGNGetSubsolver
1801
1802        """
1803        cdef TAO subsolver = TAO()
1804        CHKERR(TaoBRGNGetSubsolver(self.tao, &subsolver.tao))
1805        CHKERR(PetscINCREF(subsolver.obj))
1806        return subsolver
1807
1808    def setBRGNRegularizerObjectiveGradient(self, objgrad, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
1809        """Set the callback to compute the regularizer objective and gradient.
1810
1811        Logically collective.
1812
1813        See Also
1814        --------
1815        petsc.TaoBRGNSetRegularizerObjectiveAndGradientRoutine
1816
1817        """
1818        if args is None: args = ()
1819        if kargs is None: kargs = {}
1820        context = (objgrad, args, kargs)
1821        self.set_attr("__brgnregobjgrad__", context)
1822        CHKERR(TaoBRGNSetRegularizerObjectiveAndGradientRoutine(self.tao, TAO_BRGNRegObjGrad, <void*>context))
1823
1824    def setBRGNRegularizerHessian(self, hessian, Mat H=None, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
1825        """Set the callback to compute the regularizer Hessian.
1826
1827        Logically collective.
1828
1829        See Also
1830        --------
1831        petsc.TaoBRGNSetRegularizerHessianRoutine
1832
1833        """
1834        cdef PetscMat Hmat = NULL
1835        if H is not None: Hmat = H.mat
1836        if args is None: args = ()
1837        if kargs is None: kargs = {}
1838        context = (hessian, args, kargs)
1839        self.set_attr("__brgnreghessian__", context)
1840        CHKERR(TaoBRGNSetRegularizerHessianRoutine(self.tao, Hmat, TAO_BRGNRegHessian, <void*>context))
1841
1842    def setBRGNRegularizerWeight(self, weight: float) -> None:
1843        """Set the regularizer weight.
1844
1845        Collective.
1846
1847        """
1848        cdef PetscReal cweight = asReal(weight)
1849        CHKERR(TaoBRGNSetRegularizerWeight(self.tao, cweight))
1850
1851    def setBRGNSmoothL1Epsilon(self, epsilon: float) -> None:
1852        """Set the smooth L1 epsilon.
1853
1854        Collective.
1855
1856        See Also
1857        --------
1858        petsc.TaoBRGNSetL1SmoothEpsilon
1859
1860        """
1861        cdef PetscReal ceps = asReal(epsilon)
1862        CHKERR(TaoBRGNSetL1SmoothEpsilon(self.tao, ceps))
1863
1864    def setBRGNDictionaryMatrix(self, Mat D) -> None:
1865        """Set the dictionary matrix.
1866
1867        Collective.
1868
1869        See Also
1870        --------
1871        petsc.TaoBRGNSetDictionaryMatrix
1872
1873        """
1874        CHKERR(TaoBRGNSetDictionaryMatrix(self.tao, D.mat))
1875
1876    def getBRGNDampingVector(self) -> Vec:
1877        """Return the damping vector.
1878
1879        Not collective.
1880
1881        """
1882        # FIXME
1883        # See Also
1884        # --------
1885        # petsc.TaoBRGNGetDampingVector
1886        cdef Vec damp = Vec()
1887        CHKERR(TaoBRGNGetDampingVector(self.tao, &damp.vec))
1888        CHKERR(PetscINCREF(damp.obj))
1889        return damp
1890
1891    def createPython(self, context: Any = None, comm: Comm | None = None) -> Self:
1892        """Create an optimization solver of Python type.
1893
1894        Collective.
1895
1896        Parameters
1897        ----------
1898        context
1899            An instance of the Python class implementing the required methods.
1900        comm
1901            MPI communicator, defaults to `Sys.getDefaultComm`.
1902
1903        See Also
1904        --------
1905        petsc_python_tao, setType, setPythonContext, Type.PYTHON
1906
1907        """
1908        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_DEFAULT)
1909        cdef PetscTAO tao = NULL
1910        CHKERR(TaoCreate(ccomm, &tao))
1911        CHKERR(PetscCLEAR(self.obj)); self.tao = tao
1912        CHKERR(TaoSetType(self.tao, TAOPYTHON))
1913        CHKERR(TaoPythonSetContext(self.tao, <void*>context))
1914        return self
1915
1916    def setPythonContext(self, context: Any) -> None:
1917        """Set the instance of the class implementing the required Python methods.
1918
1919        Not collective.
1920
1921        See Also
1922        --------
1923        petsc_python_tao, getPythonContext
1924
1925        """
1926        CHKERR(TaoPythonSetContext(self.tao, <void*>context))
1927
1928    def getPythonContext(self) -> Any:
1929        """Return the instance of the class implementing the required Python methods.
1930
1931        Not collective.
1932
1933        See Also
1934        --------
1935        petsc_python_tao, setPythonContext
1936
1937        """
1938        cdef void *context = NULL
1939        CHKERR(TaoPythonGetContext(self.tao, &context))
1940        if context == NULL: return None
1941        else: return <object> context
1942
1943    def setPythonType(self, py_type: str) -> None:
1944        """Set the fully qualified Python name of the class to be used.
1945
1946        Collective.
1947
1948        See Also
1949        --------
1950        petsc_python_tao, setPythonContext, getPythonType
1951        petsc.TaoPythonSetType
1952
1953        """
1954        cdef const char *cval = NULL
1955        py_type = str2bytes(py_type, &cval)
1956        CHKERR(TaoPythonSetType(self.tao, cval))
1957
1958    def getPythonType(self) -> str:
1959        """Return the fully qualified Python name of the class used by the solver.
1960
1961        Not collective.
1962
1963        See Also
1964        --------
1965        petsc_python_tao, setPythonContext, setPythonType
1966        petsc.TaoPythonGetType
1967
1968        """
1969        cdef const char *cval = NULL
1970        CHKERR(TaoPythonGetType(self.tao, &cval))
1971        return bytes2str(cval)
1972
1973    def getLineSearch(self) -> TAOLineSearch:
1974        """Return the TAO Line Search object.
1975
1976        Not collective.
1977
1978        See Also
1979        -------
1980        petsc.TaoGetLineSearch
1981
1982        """
1983        cdef TAOLineSearch ls = TAOLineSearch()
1984        CHKERR(TaoGetLineSearch(self.tao, &ls.taols))
1985        CHKERR(PetscINCREF(ls.obj))
1986        return ls
1987
1988    # --- backward compatibility ---
1989
1990    setInitial = setSolution
1991
1992    # --- application context ---
1993
1994    property appctx:
1995        """Application context."""
1996        def __get__(self) -> Any:
1997            return self.getAppCtx()
1998
1999        def __set__(self, value: Any):
2000            self.setAppCtx(value)
2001
2002    # --- linear solver ---
2003
2004    property ksp:
2005        """Linear solver."""
2006        def __get__(self) -> KSP:
2007            return self.getKSP()
2008
2009    # --- tolerances ---
2010
2011    # FIXME: tolerances all broken
2012    property ftol:
2013        """Broken."""
2014        def __get__(self) -> Any:
2015            return self.getFunctionTolerances()
2016
2017        def __set__(self, value):
2018            if isinstance(value, (tuple, list)):
2019                self.setFunctionTolerances(*value)
2020            elif isinstance(value, dict):
2021                self.setFunctionTolerances(**value)
2022            else:
2023                raise TypeError("expecting tuple/list or dict")
2024
2025    property gtol:
2026        """Broken."""
2027        def __get__(self) -> Any:
2028            return self.getGradientTolerances()
2029
2030        def __set__(self, value):
2031            if isinstance(value, (tuple, list)):
2032                self.getGradientTolerances(*value)
2033            elif isinstance(value, dict):
2034                self.getGradientTolerances(**value)
2035            else:
2036                raise TypeError("expecting tuple/list or dict")
2037
2038    property ctol:
2039        """Broken."""
2040        def __get__(self) -> Any:
2041            return self.getConstraintTolerances()
2042
2043        def __set__(self, value):
2044            if isinstance(value, (tuple, list)):
2045                self.getConstraintTolerances(*value)
2046            elif isinstance(value, dict):
2047                self.getConstraintTolerances(**value)
2048            else:
2049                raise TypeError("expecting tuple/list or dict")
2050
2051    # --- iteration ---
2052
2053    property its:
2054        """Number of iterations."""
2055        def __get__(self) -> int:
2056            return self.getIterationNumber()
2057
2058    property gnorm:
2059        """Gradient norm."""
2060        def __get__(self) -> float:
2061            return self.getSolutionNorm()[1]
2062
2063    property cnorm:
2064        """Constraints norm."""
2065        def __get__(self) -> float:
2066            return self.getSolutionNorm()[2]
2067
2068    property solution:
2069        """Solution vector."""
2070        def __get__(self) -> Vec:
2071            return self.getSolution()
2072
2073    property objective:
2074        """Objective value."""
2075        def __get__(self) -> float:
2076            return self.getObjectiveValue()
2077
2078    property function:
2079        """Objective value."""
2080        def __get__(self) -> float:
2081            return self.getFunctionValue()
2082
2083    property gradient:
2084        """Gradient vector."""
2085        def __get__(self) -> Vec:
2086            return self.getGradient()[0]
2087
2088    # --- convergence ---
2089
2090    property reason:
2091        """Converged reason."""
2092        def __get__(self) -> ConvergedReason:
2093            return self.getConvergedReason()
2094
2095    property iterating:
2096        """Boolean indicating if the solver has not converged yet."""
2097        def __get__(self) -> bool:
2098            return self.reason == 0
2099
2100    property converged:
2101        """Boolean indicating if the solver has converged."""
2102        def __get__(self) -> bool:
2103            return self.reason > 0
2104
2105    property diverged:
2106        """Boolean indicating if the solver has failed."""
2107        def __get__(self) -> bool:
2108            return self.reason < 0
2109
2110# --------------------------------------------------------------------
2111
2112del TAOType
2113del TAOConvergedReason
2114del TAOBNCGType
2115del TAOALMMType
2116
2117# --------------------------------------------------------------------
2118
2119
2120class TAOLineSearchType:
2121    """TAO Line Search Types."""
2122    UNIT        = S_(TAOLINESEARCHUNIT)
2123    ARMIJO      = S_(TAOLINESEARCHARMIJO)
2124    MORETHUENTE = S_(TAOLINESEARCHMT)
2125    IPM         = S_(TAOLINESEARCHIPM)
2126    OWARMIJO    = S_(TAOLINESEARCHOWARMIJO)
2127    GPCG        = S_(TAOLINESEARCHGPCG)
2128
2129
2130class TAOLineSearchConvergedReason:
2131    """TAO Line Search Termination Reasons."""
2132    # iterating
2133    CONTINUE_SEARCH       = TAOLINESEARCH_CONTINUE_ITERATING
2134    # failed
2135    FAILED_INFORNAN       = TAOLINESEARCH_FAILED_INFORNAN      # inf or NaN in user function
2136    FAILED_BADPARAMETER   = TAOLINESEARCH_FAILED_BADPARAMETER  # negative value set as parameter
2137    FAILED_ASCENT         = TAOLINESEARCH_FAILED_ASCENT        # search direction is not a descent direction
2138    # succeeded
2139    SUCCESS               = TAOLINESEARCH_SUCCESS              # found step length
2140    SUCCESS_USER          = TAOLINESEARCH_SUCCESS_USER         # user-defined success criteria reached
2141    # halted
2142    HALTED_OTHER          = TAOLINESEARCH_HALTED_OTHER         # stopped search with unknown reason
2143    HALTED_MAXFCN         = TAOLINESEARCH_HALTED_MAXFCN        # maximum function evaluations reached
2144    HALTED_UPPERBOUND     = TAOLINESEARCH_HALTED_UPPERBOUND    # stopped at upper bound
2145    HALTED_LOWERBOUND     = TAOLINESEARCH_HALTED_LOWERBOUND    # stopped at lower bound
2146    HALTED_RTOL           = TAOLINESEARCH_HALTED_RTOL          # range of uncertainty is below tolerance
2147    HALTED_USER           = TAOLINESEARCH_HALTED_USER          # user-defined halt criteria reached
2148
2149# --------------------------------------------------------------------
2150
2151
2152cdef class TAOLineSearch(Object):
2153    """TAO Line Search."""
2154
2155    Type   = TAOLineSearchType
2156    ConvergedReason = TAOLineSearchConvergedReason
2157    # FIXME backward compatibility
2158    Reason = TAOLineSearchConvergedReason
2159
2160    def __cinit__(self):
2161        self.obj = <PetscObject*> &self.taols
2162        self.taols = NULL
2163
2164    def view(self, Viewer viewer=None) -> None:
2165        """View the linesearch object.
2166
2167        Collective.
2168
2169        Parameters
2170        ----------
2171        viewer
2172            A `Viewer` instance or `None` for the default viewer.
2173
2174        See Also
2175        --------
2176        petsc.TaoLineSearchView
2177
2178        """
2179        cdef PetscViewer vwr = NULL
2180        if viewer is not None: vwr = viewer.vwr
2181        CHKERR(TaoLineSearchView(self.taols, vwr))
2182
2183    def destroy(self) -> Self:
2184        """Destroy the linesearch object.
2185
2186        Collective.
2187
2188        See Also
2189        --------
2190        petsc.TaoLineSearchDestroy
2191
2192        """
2193        CHKERR(TaoLineSearchDestroy(&self.taols))
2194        return self
2195
2196    def create(self, comm=None) -> Self:
2197        """Create a TAO linesearch.
2198
2199        Collective.
2200
2201        Parameters
2202        ----------
2203        comm
2204            MPI communicator, defaults to `Sys.getDefaultComm`.
2205
2206        See Also
2207        --------
2208        Sys.getDefaultComm, petsc.TaoLineSearchCreate
2209
2210        """
2211        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_DEFAULT)
2212        cdef PetscTAOLineSearch newtaols = NULL
2213        CHKERR(TaoLineSearchCreate(ccomm, &newtaols))
2214        CHKERR(PetscCLEAR(self.obj)); self.taols = newtaols
2215        return self
2216
2217    def setType(self, ls_type: Type | str) -> None:
2218        """Set the type of the linesearch.
2219
2220        Logically collective.
2221
2222        Parameters
2223        ----------
2224        ls_type
2225            The type of the solver.
2226
2227        See Also
2228        --------
2229        getType, petsc.TaoLineSearchSetType
2230
2231        """
2232        cdef PetscTAOLineSearchType ctype = NULL
2233        ls_type = str2bytes(ls_type, &ctype)
2234        CHKERR(TaoLineSearchSetType(self.taols, ctype))
2235
2236    def getType(self) -> str:
2237        """Return the type of the linesearch.
2238
2239        Not collective.
2240
2241        See Also
2242        --------
2243        setType, petsc.TaoLineSearchGetType
2244
2245        """
2246        cdef PetscTAOLineSearchType ctype = NULL
2247        CHKERR(TaoLineSearchGetType(self.taols, &ctype))
2248        return bytes2str(ctype)
2249
2250    def setFromOptions(self) -> None:
2251        """Configure the linesearch from the options database.
2252
2253        Collective.
2254
2255        See Also
2256        --------
2257        petsc_options, petsc.TaoLineSearchSetFromOptions
2258
2259        """
2260        CHKERR(TaoLineSearchSetFromOptions(self.taols))
2261
2262    def setUp(self) -> None:
2263        """Set up the internal data structures for using the linesearch.
2264
2265        Collective.
2266
2267        See Also
2268        --------
2269        petsc.TaoLineSearchSetUp
2270
2271        """
2272        CHKERR(TaoLineSearchSetUp(self.taols))
2273
2274    def setOptionsPrefix(self, prefix: str | None = None) -> None:
2275        """Set the prefix used for searching for options in the database.
2276
2277        Logically collective.
2278
2279        See Also
2280        --------
2281        petsc_options, petsc.TaoLineSearchSetOptionsPrefix
2282
2283        """
2284        cdef const char *cprefix = NULL
2285        prefix = str2bytes(prefix, &cprefix)
2286        CHKERR(TaoLineSearchSetOptionsPrefix(self.taols, cprefix))
2287
2288    def getOptionsPrefix(self) -> str:
2289        """Return the prefix used for searching for options in the database.
2290
2291        Not collective.
2292
2293        See Also
2294        --------
2295        petsc_options, setOptionsPrefix, petsc.TaoLineSearchGetOptionsPrefix
2296
2297        """
2298        cdef const char *prefix = NULL
2299        CHKERR(TaoLineSearchGetOptionsPrefix(self.taols, &prefix))
2300        return bytes2str(prefix)
2301
2302    def setObjective(self, objective : TAOLSObjectiveFunction, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
2303        """Set the objective function evaluation callback.
2304
2305        Logically collective.
2306
2307        Parameters
2308        ----------
2309        objective
2310            The objective function callback.
2311        args
2312            Positional arguments for the callback.
2313        kargs
2314            Keyword arguments for the callback.
2315
2316        See Also
2317        --------
2318        setGradient, setObjectiveGradient
2319        petsc.TaoLineSearchSetObjectiveRoutine
2320
2321        """
2322        CHKERR(TaoLineSearchSetObjectiveRoutine(self.taols, TAOLS_Objective, NULL))
2323        if args is None: args = ()
2324        if kargs is None: kargs = {}
2325        self.set_attr("__objective__", (objective, args, kargs))
2326
2327    def setGradient(self, gradient: TAOLSGradientFunction, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
2328        """Set the gradient evaluation callback.
2329
2330        Logically collective.
2331
2332        Parameters
2333        ----------
2334        gradient
2335            The gradient callback.
2336        g
2337            The vector to store the gradient.
2338        args
2339            Positional arguments for the callback.
2340        kargs
2341            Keyword arguments for the callback.
2342
2343        See Also
2344        --------
2345        setObjective, setObjectiveGradient, setHessian
2346        petsc.TaoLineSearchSetGradientRoutine
2347
2348        """
2349        CHKERR(TaoLineSearchSetGradientRoutine(self.taols, TAOLS_Gradient, NULL))
2350        if args is None: args = ()
2351        if kargs is None: kargs = {}
2352        self.set_attr("__gradient__", (gradient, args, kargs))
2353
2354    def setObjectiveGradient(self, objgrad: TAOLSObjectiveGradientFunction, args: tuple[Any, ...] | None = None, kargs: dict[str, Any] | None = None) -> None:
2355        """Set the objective function and gradient evaluation callback.
2356
2357        Logically collective.
2358
2359        Parameters
2360        ----------
2361        objgrad
2362            The objective function and gradient callback.
2363        g
2364            The vector to store the gradient.
2365        args
2366            Positional arguments for the callback.
2367        kargs
2368            Keyword arguments for the callback.
2369
2370        See Also
2371        --------
2372        setObjective, setGradient, setHessian, getObjectiveAndGradient
2373        petsc.TaoLineSearchSetObjectiveAndGradientRoutine
2374
2375        """
2376        CHKERR(TaoLineSearchSetObjectiveAndGradientRoutine(self.taols, TAOLS_ObjGrad, NULL))
2377        if args is None: args = ()
2378        if kargs is None: kargs = {}
2379        self.set_attr("__objgrad__", (objgrad, args, kargs))
2380
2381    def useTAORoutine(self, TAO tao) -> None:
2382        """Use the objective and gradient evaluation routines from the given Tao object.
2383
2384        Logically collective.
2385
2386        See Also
2387        --------
2388        petsc.TaoLineSearchUseTaoRoutines
2389
2390        """
2391        CHKERR(TaoLineSearchUseTaoRoutines(self.taols, tao.tao))
2392
2393    def apply(self, Vec x, Vec g, Vec s) -> tuple[float, float, str]:
2394        """Performs a line-search in a given step direction.
2395
2396        Collective.
2397
2398        See Also
2399        --------
2400        petsc.TaoLineSearchApply
2401
2402        """
2403        cdef PetscReal f = 0
2404        cdef PetscReal steplen = 0
2405        cdef PetscTAOLineSearchConvergedReason reason = TAOLINESEARCH_CONTINUE_ITERATING
2406        CHKERR(TaoLineSearchApply(self.taols, x.vec, &f, g.vec, s.vec, &steplen, &reason))
2407        return (toReal(f), toReal(steplen), reason)
2408
2409    def setInitialStepLength(self, s: float) -> None:
2410        """Sets the initial step length of a line search.
2411
2412        Logically collective.
2413
2414        See Also
2415        --------
2416        petsc.TaoLineSearchSetInitialStepLength
2417
2418        """
2419        cdef PetscReal cs = asReal(s)
2420        CHKERR(TaoLineSearchSetInitialStepLength(self.taols, cs))
2421
2422# --------------------------------------------------------------------
2423
2424del TAOLineSearchType
2425del TAOLineSearchConvergedReason
2426