xref: /petsc/src/binding/petsc4py/src/petsc4py/PETSc/TAO.pyx (revision aaa8cc7d2a5c3913edcbb923e20f154fe9c4aa65)
1# --------------------------------------------------------------------
2
3class TAOType:
4    """
5    TAO Solver Types
6    """
7    LMVM     = S_(TAOLMVM)
8    NLS      = S_(TAONLS)
9    NTR      = S_(TAONTR)
10    NTL      = S_(TAONTL)
11    CG       = S_(TAOCG)
12    TRON     = S_(TAOTRON)
13    OWLQN    = S_(TAOOWLQN)
14    BMRM     = S_(TAOBMRM)
15    BLMVM    = S_(TAOBLMVM)
16    BQNLS    = S_(TAOBQNLS)
17    BNCG     = S_(TAOBNCG)
18    BNLS     = S_(TAOBNLS)
19    BNTR     = S_(TAOBNTR)
20    BNTL     = S_(TAOBNTL)
21    BQNKLS   = S_(TAOBQNKLS)
22    BQNKTR   = S_(TAOBQNKTR)
23    BQNKTL   = S_(TAOBQNKTL)
24    BQPIP    = S_(TAOBQPIP)
25    GPCG     = S_(TAOGPCG)
26    NM       = S_(TAONM)
27    POUNDERS = S_(TAOPOUNDERS)
28    BRGN     = S_(TAOBRGN)
29    LCL      = S_(TAOLCL)
30    SSILS    = S_(TAOSSILS)
31    SSFLS    = S_(TAOSSFLS)
32    ASILS    = S_(TAOASILS)
33    ASFLS    = S_(TAOASFLS)
34    IPM      = S_(TAOIPM)
35    PDIPM    = S_(TAOPDIPM)
36    SHELL    = S_(TAOSHELL)
37    ADMM     = S_(TAOADMM)
38    ALMM     = S_(TAOALMM)
39    PYTHON   = S_(TAOPYTHON)
40
41class TAOConvergedReason:
42    """
43    TAO Solver Termination Reasons
44    """
45    # iterating
46    CONTINUE_ITERATING    = TAO_CONTINUE_ITERATING    # iterating
47    CONVERGED_ITERATING   = TAO_CONTINUE_ITERATING    # iterating
48    ITERATING             = TAO_CONTINUE_ITERATING    # iterating
49    # converged
50    CONVERGED_GATOL       = TAO_CONVERGED_GATOL       # ||g(X)|| < gatol
51    CONVERGED_GRTOL       = TAO_CONVERGED_GRTOL       # ||g(X)||/f(X)  < grtol
52    CONVERGED_GTTOL       = TAO_CONVERGED_GTTOL       # ||g(X)||/||g(X0)|| < gttol
53    CONVERGED_STEPTOL     = TAO_CONVERGED_STEPTOL     # small step size
54    CONVERGED_MINF        = TAO_CONVERGED_MINF        # f(X) < F_min
55    CONVERGED_USER        = TAO_CONVERGED_USER        # user defined
56    # diverged
57    DIVERGED_MAXITS       = TAO_DIVERGED_MAXITS       #
58    DIVERGED_NAN          = TAO_DIVERGED_NAN          #
59    DIVERGED_MAXFCN       = TAO_DIVERGED_MAXFCN       #
60    DIVERGED_LS_FAILURE   = TAO_DIVERGED_LS_FAILURE   #
61    DIVERGED_TR_REDUCTION = TAO_DIVERGED_TR_REDUCTION #
62    DIVERGED_USER         = TAO_DIVERGED_USER         # user defined
63
64# --------------------------------------------------------------------
65
66cdef class TAO(Object):
67
68    """
69    TAO Solver
70    """
71
72    Type   = TAOType
73    Reason = TAOConvergedReason
74
75    def __cinit__(self):
76        self.obj = <PetscObject*> &self.tao
77        self.tao = NULL
78
79    def view(self, Viewer viewer=None):
80        """
81        """
82        cdef PetscViewer vwr = NULL
83        if viewer is not None: vwr = viewer.vwr
84        CHKERR( TaoView(self.tao, vwr) )
85
86    def destroy(self):
87        """
88        """
89        CHKERR( TaoDestroy(&self.tao) )
90        return self
91
92    def create(self, comm=None):
93        """
94        """
95        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_DEFAULT)
96        cdef PetscTAO newtao = NULL
97        CHKERR( TaoCreate(ccomm, &newtao) )
98        PetscCLEAR(self.obj); self.tao = newtao
99        return self
100
101    def setType(self, tao_type):
102        """
103        """
104        cdef PetscTAOType ctype = NULL
105        tao_type = str2bytes(tao_type, &ctype)
106        CHKERR( TaoSetType(self.tao, ctype) )
107
108    def getType(self):
109        """
110        """
111        cdef PetscTAOType ctype = NULL
112        CHKERR( TaoGetType(self.tao, &ctype) )
113        return bytes2str(ctype)
114
115    def setOptionsPrefix(self, prefix):
116        """
117        """
118        cdef const char *cprefix = NULL
119        prefix = str2bytes(prefix, &cprefix)
120        CHKERR( TaoSetOptionsPrefix(self.tao, cprefix) )
121
122    def appendOptionsPrefix(self, prefix):
123        """
124        """
125        cdef const char *cprefix = NULL
126        prefix = str2bytes(prefix, &cprefix)
127        CHKERR( TaoAppendOptionsPrefix(self.tao, cprefix) )
128
129    def getOptionsPrefix(self):
130        """
131        """
132        cdef const char *prefix = NULL
133        CHKERR( TaoGetOptionsPrefix(self.tao, &prefix) )
134        return bytes2str(prefix)
135
136    def setFromOptions(self):
137        """
138        """
139        CHKERR( TaoSetFromOptions(self.tao) )
140
141    def setUp(self):
142        """
143        """
144        CHKERR( TaoSetUp(self.tao) )
145
146    #
147
148    def setInitialTrustRegionRadius(self, radius):
149        cdef PetscReal cradius = asReal(radius)
150        CHKERR( TaoSetInitialTrustRegionRadius(self.tao, cradius) )
151
152    # --------------
153
154    def setAppCtx(self, appctx):
155        self.set_attr("__appctx__", appctx)
156
157    def getAppCtx(self):
158        return self.get_attr("__appctx__")
159
160    def setSolution(self, Vec x):
161        """
162        """
163        CHKERR( TaoSetSolution(self.tao, x.vec) )
164
165    def setObjective(self, objective, args=None, kargs=None):
166        """
167        """
168        if args is None: args = ()
169        if kargs is None: kargs = {}
170        context = (objective, args, kargs)
171        self.set_attr("__objective__", context)
172        CHKERR( TaoSetObjective(self.tao, TAO_Objective, <void*>context) )
173
174    def setResidual(self, residual, Vec R=None, args=None, kargs=None):
175        """
176        """
177        cdef PetscVec Rvec = NULL
178        if R is not None: Rvec = R.vec
179        if args is None: args = ()
180        if kargs is None: kargs = {}
181        context = (residual, args, kargs)
182        self.set_attr("__residual__", context)
183        CHKERR( TaoSetResidualRoutine(self.tao, Rvec, TAO_Residual, <void*>context) )
184
185    def setJacobianResidual(self, jacobian, Mat J=None, Mat P=None, args=None, kargs=None):
186        """
187        """
188        cdef PetscMat Jmat = NULL
189        if J is not None: Jmat = J.mat
190        cdef PetscMat Pmat = Jmat
191        if P is not None: Pmat = P.mat
192        if args is None: args = ()
193        if kargs is None: kargs = {}
194        context = (jacobian, args, kargs)
195        self.set_attr("__jacobian_residual__", context)
196        CHKERR( TaoSetJacobianResidualRoutine(self.tao, Jmat, Pmat, TAO_JacobianResidual, <void*>context) )
197
198    def setGradient(self, gradient, Vec g or None, args=None, kargs=None):
199        """
200        """
201        cdef PetscVec gvec = NULL
202        if g is not None: gvec = g.vec
203        if args is None: args = ()
204        if kargs is None: kargs = {}
205        context = (gradient, args, kargs)
206        self.set_attr("__gradient__", context)
207        CHKERR( TaoSetGradient(self.tao, gvec, TAO_Gradient, <void*>context) )
208
209    def getGradient(self):
210        """
211        """
212        cdef Vec vec = Vec()
213        CHKERR( TaoGetGradient(self.tao, &vec.vec, NULL, NULL) )
214        PetscINCREF(vec.obj)
215        cdef object gradient = self.get_attr("__gradient__")
216        return (vec, gradient)
217
218    def setObjectiveGradient(self, objgrad, Vec g or None, args=None, kargs=None):
219        """
220        """
221        cdef PetscVec gvec = NULL
222        if g is not None: gvec = g.vec
223        if args is None: args = ()
224        if kargs is None: kargs = {}
225        context = (objgrad, args, kargs)
226        self.set_attr("__objgrad__", context)
227        CHKERR( TaoSetObjectiveAndGradient(self.tao, gvec, TAO_ObjGrad, <void*>context) )
228
229    def getObjectiveAndGradient(self):
230        """
231        """
232        cdef Vec vec = Vec()
233        CHKERR( TaoGetObjectiveAndGradient(self.tao, &vec.vec, NULL, NULL) )
234        PetscINCREF(vec.obj)
235        cdef object objgrad = self.get_attr("__objgrad__")
236        return (vec, objgrad)
237
238    def setVariableBounds(self, varbounds, args=None, kargs=None):
239        """
240        """
241        cdef Vec xl = None, xu = None
242        if (isinstance(varbounds, list) or isinstance(varbounds, tuple)):
243            ol, ou = varbounds
244            xl = <Vec?> ol; xu = <Vec?> ou
245            CHKERR( TaoSetVariableBounds(self.tao, xl.vec, xu.vec) )
246            return
247        if isinstance(varbounds, Vec):
248            ol = varbounds; ou = args
249            xl = <Vec?> ol; xu = <Vec?> ou
250            CHKERR( TaoSetVariableBounds(self.tao, xl.vec, xu.vec) )
251            return
252        if args is None: args = ()
253        if kargs is None: kargs = {}
254        context = (varbounds, args, kargs)
255        self.set_attr("__varbounds__", context)
256        CHKERR( TaoSetVariableBoundsRoutine(self.tao, TAO_VarBounds, <void*>context) )
257
258    def setConstraints(self, constraints, Vec C=None, args=None, kargs=None):
259        """
260        """
261        cdef PetscVec Cvec = NULL
262        if C is not None: Cvec = C.vec
263        if args is None: args = ()
264        if kargs is None: kargs = {}
265        context = (constraints, args, kargs)
266        self.set_attr("__constraints__", context)
267        CHKERR( TaoSetConstraintsRoutine(self.tao, Cvec, TAO_Constraints, <void*>context) )
268
269    def setHessian(self, hessian, Mat H=None, Mat P=None,
270                   args=None, kargs=None):
271        cdef PetscMat Hmat = NULL
272        if H is not None: Hmat = H.mat
273        cdef PetscMat Pmat = Hmat
274        if P is not None: Pmat = P.mat
275        if args is None: args = ()
276        if kargs is None: kargs = {}
277        context = (hessian, args, kargs)
278        self.set_attr("__hessian__", context)
279        CHKERR( TaoSetHessian(self.tao, Hmat, Pmat, TAO_Hessian, <void*>context) )
280
281    def getHessian(self):
282        cdef Mat J = Mat()
283        cdef Mat P = Mat()
284        CHKERR( TaoGetHessian(self.tao, &J.mat, &P.mat, NULL, NULL) )
285        PetscINCREF(J.obj)
286        PetscINCREF(P.obj)
287        cdef object hessian = self.get_attr("__hessian__")
288        return (J, P, hessian)
289
290    def setJacobian(self, jacobian, Mat J=None, Mat P=None,
291                    args=None, kargs=None):
292        """
293        """
294        cdef PetscMat Jmat = NULL
295        if J is not None: Jmat = J.mat
296        cdef PetscMat Pmat = Jmat
297        if P is not None: Pmat = P.mat
298        if args is None: args = ()
299        if kargs is None: kargs = {}
300        context = (jacobian, args, kargs)
301        self.set_attr("__jacobian__", context)
302        CHKERR( TaoSetJacobianRoutine(self.tao, Jmat, Pmat, TAO_Jacobian, <void*>context) )
303
304    #
305
306    def setStateDesignIS(self, IS state=None, IS design=None):
307        """
308        """
309        cdef PetscIS s_is = NULL, d_is = NULL
310        if state  is not None: s_is = state.iset
311        if design is not None: d_is = design.iset
312        CHKERR( TaoSetStateDesignIS(self.tao, s_is, d_is) )
313
314    def setJacobianState(self, jacobian_state, Mat J=None, Mat P=None, Mat I=None,
315                         args=None, kargs=None):
316        """
317        """
318        cdef PetscMat Jmat = NULL
319        if J is not None: Jmat = J.mat
320        cdef PetscMat Pmat = Jmat
321        if P is not None: Pmat = P.mat
322        cdef PetscMat Imat = NULL
323        if I is not None: Imat = I.mat
324        if args is None: args = ()
325        if kargs is None: kargs = {}
326        context = (jacobian_state, args, kargs)
327        self.set_attr("__jacobian_state__", context)
328        CHKERR( TaoSetJacobianStateRoutine(self.tao, Jmat, Pmat, Imat,
329                                           TAO_JacobianState, <void*>context) )
330
331    def setJacobianDesign(self, jacobian_design, Mat J=None,
332                          args=None, kargs=None):
333        """
334        """
335        cdef PetscMat Jmat = NULL
336        if J is not None: Jmat = J.mat
337        if args is None: args = ()
338        if kargs is None: kargs = {}
339        context = (jacobian_design, args, kargs)
340        self.set_attr("__jacobian_design__", context)
341        CHKERR( TaoSetJacobianDesignRoutine(self.tao, Jmat,
342                                            TAO_JacobianDesign, <void*>context) )
343
344
345    def setEqualityConstraints(self, equality_constraints, Vec c,
346                               args=None, kargs=None):
347        """
348        """
349        if args is None: args = ()
350        if kargs is None: kargs = {}
351        context = (equality_constraints, args, kargs)
352        self.set_attr("__equality_constraints__", context)
353        CHKERR( TaoSetEqualityConstraintsRoutine(self.tao, c.vec,
354                                                 TAO_EqualityConstraints, <void*>context) )
355
356
357    def setJacobianEquality(self, jacobian_equality, Mat J=None, Mat P=None,
358                            args=None, kargs=None):
359        """
360        """
361        cdef PetscMat Jmat = NULL
362        if J is not None: Jmat = J.mat
363        cdef PetscMat Pmat = Jmat
364        if P is not None: Pmat = P.mat
365        if args is None: args = ()
366        if kargs is None: kargs = {}
367        context = (jacobian_equality, args, kargs)
368        self.set_attr("__jacobian_equality__", context)
369        CHKERR( TaoSetJacobianEqualityRoutine(self.tao, Jmat, Pmat,
370                                              TAO_JacobianEquality, <void*>context) )
371
372    def setUpdate(self, update, args=None, kargs=None):
373        """
374        """
375        if update is not None:
376            if args  is None: args  = ()
377            if kargs is None: kargs = {}
378            context = (update, args, kargs)
379            self.set_attr('__update__', context)
380            CHKERR( TaoSetUpdate(self.tao, TAO_Update, <void*>context) )
381        else:
382            self.set_attr('__update__', None)
383            CHKERR( TaoSetUpdate(self.tao, NULL, NULL) )
384
385    def getUpdate(self):
386        return self.get_attr('__update__')
387
388    # --------------
389
390    def computeObjective(self, Vec x):
391        """
392        """
393        cdef PetscReal f = 0
394        CHKERR( TaoComputeObjective(self.tao, x.vec, &f) )
395        return toReal(f)
396
397    def computeResidual(self, Vec x, Vec f):
398        """
399        """
400        CHKERR( TaoComputeResidual(self.tao, x.vec, f.vec) )
401
402    def computeGradient(self, Vec x, Vec g):
403        """
404        """
405        CHKERR( TaoComputeGradient(self.tao, x.vec, g.vec) )
406
407    def computeObjectiveGradient(self, Vec x, Vec g):
408        """
409        """
410        cdef PetscReal f = 0
411        CHKERR( TaoComputeObjectiveAndGradient(self.tao, x.vec, &f, g.vec) )
412        return toReal(f)
413
414    def computeDualVariables(self, Vec xl, Vec xu):
415        """
416        """
417        CHKERR( TaoComputeDualVariables(self.tao, xl.vec, xu.vec) )
418
419    def computeVariableBounds(self, Vec xl, Vec xu):
420        """
421        """
422        CHKERR( TaoComputeVariableBounds(self.tao) )
423        cdef PetscVec Lvec = NULL, Uvec = NULL
424        CHKERR( TaoGetVariableBounds(self.tao, &Lvec, &Uvec) )
425        if xl.vec != NULL:
426            if Lvec != NULL:
427                CHKERR( VecCopy(Lvec, xl.vec) )
428            else:
429                CHKERR( VecSet(xl.vec, <PetscScalar>PETSC_NINFINITY) )
430        if xu.vec != NULL:
431            if Uvec != NULL:
432                CHKERR( VecCopy(Uvec, xu.vec) )
433            else:
434                CHKERR( VecSet(xu.vec, <PetscScalar>PETSC_INFINITY) )
435
436    def computeConstraints(self, Vec x, Vec c):
437        """
438        """
439        CHKERR( TaoComputeConstraints(self.tao, x.vec, c.vec) )
440
441    def computeHessian(self, Vec x, Mat H, Mat P=None):
442        """
443        """
444        cdef PetscMat hmat = H.mat, pmat = H.mat
445        if P is not None: pmat = P.mat
446        CHKERR( TaoComputeHessian(self.tao, x.vec, hmat, pmat) )
447
448    def computeJacobian(self, Vec x, Mat J, Mat P=None):
449        """
450        """
451        cdef PetscMat jmat = J.mat, pmat = J.mat
452        if P is not None: pmat = P.mat
453        CHKERR( TaoComputeJacobian(self.tao, x.vec, jmat, pmat) )
454
455    # --------------
456
457    #
458
459    def setTolerances(self, gatol=None, grtol=None, gttol=None):
460        """
461        """
462        cdef PetscReal _gatol=PETSC_DEFAULT, _grtol=PETSC_DEFAULT, _gttol=PETSC_DEFAULT
463        if gatol is not None: _gatol = asReal(gatol)
464        if grtol is not None: _grtol = asReal(grtol)
465        if gttol is not None: _gttol = asReal(gttol)
466        CHKERR( TaoSetTolerances(self.tao, _gatol, _grtol, _gttol) )
467
468    def getTolerances(self):
469        """
470        """
471        cdef PetscReal _gatol=PETSC_DEFAULT, _grtol=PETSC_DEFAULT, _gttol=PETSC_DEFAULT
472        CHKERR( TaoGetTolerances(self.tao, &_gatol, &_grtol, &_gttol) )
473        return (toReal(_gatol), toReal(_grtol), toReal(_gttol))
474
475    def setMaximumIterations(self, mit):
476        """
477        """
478        cdef PetscInt _mit = asInt(mit)
479        CHKERR( TaoSetMaximumIterations(self.tao, _mit) )
480
481    def getMaximumIterations(self):
482        """
483        """
484        cdef PetscInt _mit = PETSC_DEFAULT
485        CHKERR( TaoGetMaximumIterations(self.tao, &_mit) )
486        return toInt(_mit)
487
488    def setMaximumFunctionEvaluations(self, mit):
489        """
490        """
491        cdef PetscInt _mit = asInt(mit)
492        CHKERR( TaoSetMaximumFunctionEvaluations(self.tao, _mit) )
493
494    def getMaximumFunctionEvaluations(self):
495        """
496        """
497        cdef PetscInt _mit = PETSC_DEFAULT
498        CHKERR( TaoGetMaximumFunctionEvaluations(self.tao, &_mit) )
499        return toInt(_mit)
500
501    def setConstraintTolerances(self, catol=None, crtol=None):
502        """
503        """
504        cdef PetscReal _catol=PETSC_DEFAULT, _crtol=PETSC_DEFAULT
505        if catol is not None: _catol = asReal(catol)
506        if crtol is not None: _crtol = asReal(crtol)
507        CHKERR( TaoSetConstraintTolerances(self.tao, _catol, _crtol) )
508
509    def getConstraintTolerances(self):
510        """
511        """
512        cdef PetscReal _catol=PETSC_DEFAULT, _crtol=PETSC_DEFAULT
513        CHKERR( TaoGetConstraintTolerances(self.tao, &_catol, &_crtol) )
514        return (toReal(_catol), toReal(_crtol))
515
516    def setConvergenceTest(self, converged, args=None, kargs=None):
517        """
518        """
519        if converged is None:
520            CHKERR( TaoSetConvergenceTest(self.tao, TaoDefaultConvergenceTest, NULL) )
521            self.set_attr('__converged__', None)
522        else:
523            if args is None: args = ()
524            if kargs is None: kargs = {}
525            self.set_attr('__converged__', (converged, args, kargs))
526            CHKERR( TaoSetConvergenceTest(self.tao, TAO_Converged, NULL) )
527
528    def getConvergenceTest(self):
529        """
530        """
531        return self.get_attr('__converged__')
532
533    def setConvergedReason(self, reason):
534        """
535        """
536        cdef PetscTAOConvergedReason creason = reason
537        CHKERR( TaoSetConvergedReason(self.tao, creason) )
538
539    def getConvergedReason(self):
540        """
541        """
542        cdef PetscTAOConvergedReason creason = TAO_CONTINUE_ITERATING
543        CHKERR( TaoGetConvergedReason(self.tao, &creason) )
544        return creason
545
546    def setMonitor(self, monitor, args=None, kargs=None):
547        """
548        """
549        if monitor is None: return
550        cdef object monitorlist = self.get_attr('__monitor__')
551        if args  is None: args  = ()
552        if kargs is None: kargs = {}
553        if monitorlist is None:
554            CHKERR( TaoSetMonitor(self.tao, TAO_Monitor, NULL, NULL) )
555            self.set_attr('__monitor__',  [(monitor, args, kargs)])
556        else:
557            monitorlist.append((monitor, args, kargs))
558
559    def getMonitor(self):
560        """
561        """
562        return self.get_attr('__monitor__')
563
564    def cancelMonitor(self):
565        """
566        """
567        CHKERR( TaoCancelMonitors(self.tao) )
568        self.set_attr('__monitor__',  None)
569
570    # Tao overwrites this statistics. Copy user defined only if present
571    def monitor(self, its=None, f=None, res=None, cnorm=None, step=None):
572        cdef PetscInt cits = 0
573        cdef PetscReal cf = 0.0
574        cdef PetscReal cres = 0.0
575        cdef PetscReal ccnorm = 0.0
576        cdef PetscReal cstep = 0.0
577        CHKERR( TaoGetSolutionStatus(self.tao, &cits, &cf, &cres, &ccnorm, &cstep, NULL) )
578        if its is not None:
579            cits = asInt(its)
580        if f is not None:
581            cf = asReal(f)
582        if res is not None:
583            cres = asReal(res)
584        if cnorm is not None:
585            ccnorm = asReal(cnorm)
586        if step is not None:
587            cstep = asReal(step)
588        CHKERR( TaoMonitor(self.tao, cits, cf, cres, ccnorm, cstep) )
589
590    #
591
592    def solve(self, Vec x=None):
593        """
594        """
595        if x is not None:
596            CHKERR( TaoSetSolution(self.tao, x.vec) )
597        CHKERR( TaoSolve(self.tao) )
598
599    def getSolution(self):
600        """
601        """
602        cdef Vec vec = Vec()
603        CHKERR( TaoGetSolution(self.tao, &vec.vec) )
604        PetscINCREF(vec.obj)
605        return vec
606
607    def setGradientNorm(self, Mat mat):
608        """
609        """
610        CHKERR( TaoSetGradientNorm(self.tao, mat.mat) )
611
612    def getGradientNorm(self):
613        """
614        """
615        cdef Mat mat = Mat()
616        CHKERR( TaoGetGradientNorm(self.tao, &mat.mat) )
617        PetscINCREF(mat.obj)
618        return mat
619
620    def setLMVMH0(self, Mat mat):
621        """
622        """
623        CHKERR( TaoLMVMSetH0(self.tao, mat.mat) )
624
625    def getLMVMH0(self):
626        """
627        """
628        cdef Mat mat = Mat()
629        CHKERR( TaoLMVMGetH0(self.tao, &mat.mat) )
630        PetscINCREF(mat.obj)
631        return mat
632
633    def getLMVMH0KSP(self):
634        """
635        """
636        cdef KSP ksp = KSP()
637        CHKERR( TaoLMVMGetH0KSP(self.tao, &ksp.ksp) )
638        PetscINCREF(ksp.obj)
639        return ksp
640
641    def getVariableBounds(self):
642        """
643        """
644        cdef Vec xl = Vec(), xu = Vec()
645        CHKERR( TaoGetVariableBounds(self.tao, &xl.vec, &xu.vec) )
646        PetscINCREF(xl.obj); PetscINCREF(xu.obj)
647        return (xl, xu)
648
649    def setIterationNumber(self, its):
650        """
651        """
652        cdef PetscInt ival = asInt(its)
653        CHKERR( TaoSetIterationNumber(self.tao, ival) )
654
655    def getIterationNumber(self):
656        """
657        """
658        cdef PetscInt its=0
659        CHKERR( TaoGetIterationNumber(self.tao, &its) )
660        return toInt(its)
661
662    def getObjectiveValue(self):
663        """
664        """
665        cdef PetscReal fval=0
666        CHKERR( TaoGetSolutionStatus(self.tao, NULL, &fval, NULL, NULL, NULL, NULL) )
667        return toReal(fval)
668
669    getFunctionValue = getObjectiveValue
670
671    def getConvergedReason(self):
672        """
673        """
674        cdef PetscTAOConvergedReason reason = TAO_CONTINUE_ITERATING
675        CHKERR( TaoGetConvergedReason(self.tao, &reason) )
676        return reason
677
678    def getSolutionNorm(self):
679        """
680        """
681        cdef PetscReal gnorm=0
682        cdef PetscReal cnorm=0
683        cdef PetscReal fval=0
684        CHKERR( TaoGetSolutionStatus(self.tao, NULL, &fval, &gnorm, &cnorm, NULL, NULL) )
685        return (toReal(fval), toReal(gnorm), toReal(cnorm))
686
687    def getSolutionStatus(self):
688        """
689        """
690        cdef PetscInt its=0
691        cdef PetscReal fval=0, gnorm=0, cnorm=0, xdiff=0
692        cdef PetscTAOConvergedReason reason = TAO_CONTINUE_ITERATING
693        CHKERR( TaoGetSolutionStatus(self.tao, &its,
694                                     &fval, &gnorm, &cnorm, &xdiff,
695                                     &reason) )
696        return (toInt(its), toReal(fval),
697                toReal(gnorm), toReal(cnorm),
698                toReal(xdiff), reason)
699
700    def getKSP(self):
701        """
702        """
703        cdef KSP ksp = KSP()
704        CHKERR( TaoGetKSP(self.tao, &ksp.ksp) )
705        PetscINCREF(ksp.obj)
706        return ksp
707
708    # BRGN routines
709
710    def getBRGNSubsolver(self):
711        """
712        """
713        cdef TAO subsolver = TAO()
714        CHKERR( TaoBRGNGetSubsolver(self.tao, &subsolver.tao) )
715        PetscINCREF(subsolver.obj)
716        return subsolver
717
718    def setBRGNRegularizerObjectiveGradient(self, objgrad, args=None, kargs=None):
719        """
720        """
721        if args is None: args = ()
722        if kargs is None: kargs = {}
723        context = (objgrad, args, kargs)
724        self.set_attr("__brgnregobjgrad__", context)
725        CHKERR( TaoBRGNSetRegularizerObjectiveAndGradientRoutine(self.tao, TAO_BRGNRegObjGrad, <void*>context) )
726
727    def setBRGNRegularizerHessian(self, hessian, Mat H=None, args=None, kargs=None):
728        cdef PetscMat Hmat = NULL
729        if H is not None: Hmat = H.mat
730        if args is None: args = ()
731        if kargs is None: kargs = {}
732        context = (hessian, args, kargs)
733        self.set_attr("__brgnreghessian__", context)
734        CHKERR( TaoBRGNSetRegularizerHessianRoutine(self.tao, Hmat, TAO_BRGNRegHessian, <void*>context) )
735
736    def setBRGNRegularizerWeight(self, weight):
737        """
738        """
739        cdef PetscReal cweight = asReal(weight)
740        CHKERR( TaoBRGNSetRegularizerWeight(self.tao, cweight) )
741
742    def setBRGNSmoothL1Epsilon(self, epsilon):
743        """
744        """
745        cdef PetscReal ceps = asReal(epsilon)
746        CHKERR( TaoBRGNSetL1SmoothEpsilon(self.tao, ceps) )
747
748    def setBRGNDictionaryMatrix(self, Mat D):
749        """
750        """
751        CHKERR( TaoBRGNSetDictionaryMatrix(self.tao, D.mat) )
752
753    def getBRGNDampingVector(self):
754        """
755        """
756        cdef Vec damp = Vec()
757        CHKERR( TaoBRGNGetDampingVector(self.tao, &damp.vec) )
758        PetscINCREF(damp.obj)
759        return damp
760
761    def createPython(self, context=None, comm=None):
762        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_DEFAULT)
763        cdef PetscTAO tao = NULL
764        CHKERR( TaoCreate(ccomm, &tao) )
765        PetscCLEAR(self.obj); self.tao = tao
766        CHKERR( TaoSetType(self.tao, TAOPYTHON) )
767        CHKERR( TaoPythonSetContext(self.tao, <void*>context) )
768        return self
769
770    def setPythonContext(self, context):
771        CHKERR( TaoPythonSetContext(self.tao, <void*>context) )
772
773    def getPythonContext(self):
774        cdef void *context = NULL
775        CHKERR( TaoPythonGetContext(self.tao, &context) )
776        if context == NULL: return None
777        else: return <object> context
778
779    def setPythonType(self, py_type):
780        cdef const char *cval = NULL
781        py_type = str2bytes(py_type, &cval)
782        CHKERR( TaoPythonSetType(self.tao, cval) )
783
784    def getPythonType(self):
785        cdef const char *cval = NULL
786        CHKERR( TaoPythonGetType(self.tao, &cval) )
787        return bytes2str(cval)
788
789    # --- backward compatibility ---
790
791    setInitial = setSolution
792
793    # --- application context ---
794
795    property appctx:
796        def __get__(self):
797            return self.getAppCtx()
798        def __set__(self, value):
799            self.setAppCtx(value)
800
801    # --- linear solver ---
802
803    property ksp:
804        def __get__(self):
805            return self.getKSP()
806
807    # --- tolerances ---
808
809    property ftol:
810        def __get__(self):
811            return self.getFunctionTolerances()
812        def __set__(self, value):
813            if isinstance(value, (tuple, list)):
814                self.setFunctionTolerances(*value)
815            elif isinstance(value, dict):
816                self.setFunctionTolerances(**value)
817            else:
818                raise TypeError("expecting tuple/list or dict")
819
820    property gtol:
821        def __get__(self):
822            return self.getGradientTolerances()
823        def __set__(self, value):
824            if isinstance(value, (tuple, list)):
825                self.getGradientTolerances(*value)
826            elif isinstance(value, dict):
827                self.getGradientTolerances(**value)
828            else:
829                raise TypeError("expecting tuple/list or dict")
830
831    property ctol:
832        def __get__(self):
833            return self.getConstraintTolerances()
834        def __set__(self, value):
835            if isinstance(value, (tuple, list)):
836                self.getConstraintTolerances(*value)
837            elif isinstance(value, dict):
838                self.getConstraintTolerances(**value)
839            else:
840                raise TypeError("expecting tuple/list or dict")
841
842    # --- iteration ---
843
844    property its:
845        def __get__(self):
846            return self.getIterationNumber()
847
848    property gnorm:
849        def __get__(self):
850            return self.getSolutionNorm()[1]
851
852    property cnorm:
853        def __get__(self):
854            return self.getSolutionNorm()[2]
855
856    property solution:
857        def __get__(self):
858            return self.getSolution()
859
860    property objective:
861        def __get__(self):
862            return self.getObjectiveValue()
863
864    property function:
865        def __get__(self):
866            return self.getFunctionValue()
867
868    property gradient:
869        def __get__(self):
870            return self.getGradient()[0]
871
872    # --- convergence ---
873
874    property reason:
875        def __get__(self):
876            return self.getConvergedReason()
877
878    property iterating:
879        def __get__(self):
880            return self.reason == 0
881
882    property converged:
883        def __get__(self):
884            return self.reason > 0
885
886    property diverged:
887        def __get__(self):
888            return self.reason < 0
889
890# --------------------------------------------------------------------
891
892del TAOType
893del TAOConvergedReason
894
895# --------------------------------------------------------------------
896