xref: /petsc/src/tao/interface/taosolver.c (revision 40cbb1a031ea8f2be4fe2b92dc842b003ad37be3)
1 #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/
2 #include <petsc/private/snesimpl.h>
3 
4 PetscBool TaoRegisterAllCalled = PETSC_FALSE;
5 PetscFunctionList TaoList = NULL;
6 
7 PetscClassId TAO_CLASSID;
8 
9 PetscLogEvent TAO_Solve;
10 PetscLogEvent TAO_ObjectiveEval;
11 PetscLogEvent TAO_GradientEval;
12 PetscLogEvent TAO_ObjGradEval;
13 PetscLogEvent TAO_HessianEval;
14 PetscLogEvent TAO_JacobianEval;
15 PetscLogEvent TAO_ConstraintsEval;
16 
17 const char *TaoSubSetTypes[] = {"subvec","mask","matrixfree","TaoSubSetType","TAO_SUBSET_",NULL};
18 
19 struct _n_TaoMonitorDrawCtx {
20   PetscViewer viewer;
21   PetscInt    howoften;  /* when > 0 uses iteration % howoften, when negative only final solution plotted */
22 };
23 
24 static PetscErrorCode KSPPreSolve_TAOEW_Private(KSP ksp, Vec b, Vec x, Tao tao)
25 {
26   SNES snes_ewdummy = tao->snes_ewdummy;
27 
28   PetscFunctionBegin;
29   if (!snes_ewdummy) PetscFunctionReturn(0);
30   /* populate snes_ewdummy struct values used in KSPPreSolve_SNESEW */
31   snes_ewdummy->vec_func = b;
32   snes_ewdummy->rtol = tao->gttol;
33   snes_ewdummy->iter = tao->niter;
34   PetscCall(VecNorm(b,NORM_2,&snes_ewdummy->norm));
35   PetscCall(KSPPreSolve_SNESEW(ksp,b,x,snes_ewdummy));
36   snes_ewdummy->vec_func = NULL;
37   PetscFunctionReturn(0);
38 }
39 
40 static PetscErrorCode KSPPostSolve_TAOEW_Private(KSP ksp, Vec b, Vec x, Tao tao)
41 {
42   SNES snes_ewdummy = tao->snes_ewdummy;
43 
44   PetscFunctionBegin;
45   if (!snes_ewdummy) PetscFunctionReturn(0);
46   PetscCall(KSPPostSolve_SNESEW(ksp,b,x,snes_ewdummy));
47   PetscFunctionReturn(0);
48 }
49 
50 static PetscErrorCode TaoSetUpEW_Private(Tao tao)
51 {
52   SNESKSPEW  *kctx;
53   const char *ewprefix;
54 
55   PetscFunctionBegin;
56   if (!tao->ksp) PetscFunctionReturn(0);
57   if (tao->ksp_ewconv) {
58     if (!tao->snes_ewdummy) PetscCall(SNESCreate(PetscObjectComm((PetscObject)tao),&tao->snes_ewdummy));
59     tao->snes_ewdummy->ksp_ewconv = PETSC_TRUE;
60     PetscCall(KSPSetPreSolve(tao->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPreSolve_TAOEW_Private,tao));
61     PetscCall(KSPSetPostSolve(tao->ksp,(PetscErrorCode (*)(KSP,Vec,Vec,void*))KSPPostSolve_TAOEW_Private,tao));
62 
63     PetscCall(KSPGetOptionsPrefix(tao->ksp,&ewprefix));
64     kctx = (SNESKSPEW*)tao->snes_ewdummy->kspconvctx;
65     PetscCall(SNESEWSetFromOptions_Private(kctx,PetscObjectComm((PetscObject)tao),ewprefix));
66   } else PetscCall(SNESDestroy(&tao->snes_ewdummy));
67   PetscFunctionReturn(0);
68 }
69 
70 /*@
71   TaoCreate - Creates a TAO solver
72 
73   Collective
74 
75   Input Parameter:
76 . comm - MPI communicator
77 
78   Output Parameter:
79 . newtao - the new Tao context
80 
81   Available methods include:
82 +    nls - Newton's method with line search for unconstrained minimization
83 .    ntr - Newton's method with trust region for unconstrained minimization
84 .    ntl - Newton's method with trust region, line search for unconstrained minimization
85 .    lmvm - Limited memory variable metric method for unconstrained minimization
86 .    cg - Nonlinear conjugate gradient method for unconstrained minimization
87 .    nm - Nelder-Mead algorithm for derivate-free unconstrained minimization
88 .    tron - Newton Trust Region method for bound constrained minimization
89 .    gpcg - Newton Trust Region method for quadratic bound constrained minimization
90 .    blmvm - Limited memory variable metric method for bound constrained minimization
91 .    lcl - Linearly constrained Lagrangian method for pde-constrained minimization
92 -    pounders - Model-based algorithm for nonlinear least squares
93 
94    Options Database Keys:
95 .   -tao_type - select which method TAO should use
96 
97    Level: beginner
98 
99 .seealso: `TaoSolve()`, `TaoDestroy()`
100 @*/
101 PetscErrorCode TaoCreate(MPI_Comm comm, Tao *newtao)
102 {
103   Tao            tao;
104 
105   PetscFunctionBegin;
106   PetscValidPointer(newtao,2);
107   PetscCall(TaoInitializePackage());
108   PetscCall(TaoLineSearchInitializePackage());
109   PetscCall(PetscHeaderCreate(tao,TAO_CLASSID,"Tao","Optimization solver","Tao",comm,TaoDestroy,TaoView));
110 
111   /* Set non-NULL defaults */
112   tao->ops->convergencetest = TaoDefaultConvergenceTest;
113 
114   tao->max_it    = 10000;
115   tao->max_funcs = -1;
116 #if defined(PETSC_USE_REAL_SINGLE)
117   tao->gatol     = 1e-5;
118   tao->grtol     = 1e-5;
119   tao->crtol     = 1e-5;
120   tao->catol     = 1e-5;
121 #else
122   tao->gatol     = 1e-8;
123   tao->grtol     = 1e-8;
124   tao->crtol     = 1e-8;
125   tao->catol     = 1e-8;
126 #endif
127   tao->gttol     = 0.0;
128   tao->steptol   = 0.0;
129   tao->trust0    = PETSC_INFINITY;
130   tao->fmin      = PETSC_NINFINITY;
131 
132   tao->hist_reset  = PETSC_TRUE;
133 
134   PetscCall(TaoResetStatistics(tao));
135   *newtao = tao;
136   PetscFunctionReturn(0);
137 }
138 
139 /*@
140   TaoSolve - Solves an optimization problem min F(x) s.t. l <= x <= u
141 
142   Collective on Tao
143 
144   Input Parameters:
145 . tao - the Tao context
146 
147   Notes:
148   The user must set up the Tao with calls to TaoSetSolution(),
149   TaoSetObjective(),
150   TaoSetGradient(), and (if using 2nd order method) TaoSetHessian().
151 
152   You should call TaoGetConvergedReason() or run with -tao_converged_reason to determine if the optimization algorithm actually succeeded or
153   why it failed.
154 
155   Level: beginner
156 
157 .seealso: `TaoCreate()`, `TaoSetObjective()`, `TaoSetGradient()`, `TaoSetHessian()`, `TaoGetConvergedReason()`
158  @*/
159 PetscErrorCode TaoSolve(Tao tao)
160 {
161   static PetscBool set = PETSC_FALSE;
162 
163   PetscFunctionBegin;
164   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
165   PetscCall(PetscCitationsRegister("@TechReport{tao-user-ref,\n"
166                                    "title   = {Toolkit for Advanced Optimization (TAO) Users Manual},\n"
167                                    "author  = {Todd Munson and Jason Sarich and Stefan Wild and Steve Benson and Lois Curfman McInnes},\n"
168                                    "Institution = {Argonne National Laboratory},\n"
169                                    "Year   = 2014,\n"
170                                    "Number = {ANL/MCS-TM-322 - Revision 3.5},\n"
171                                    "url    = {https://www.mcs.anl.gov/research/projects/tao/}\n}\n",&set));
172   tao->header_printed = PETSC_FALSE;
173   PetscCall(TaoSetUp(tao));
174   PetscCall(TaoResetStatistics(tao));
175   if (tao->linesearch) {
176     PetscCall(TaoLineSearchReset(tao->linesearch));
177   }
178 
179   PetscCall(PetscLogEventBegin(TAO_Solve,tao,0,0,0));
180   if (tao->ops->solve) PetscCall((*tao->ops->solve)(tao));
181   PetscCall(PetscLogEventEnd(TAO_Solve,tao,0,0,0));
182 
183   PetscCall(VecViewFromOptions(tao->solution,(PetscObject)tao,"-tao_view_solution"));
184 
185   tao->ntotalits += tao->niter;
186   PetscCall(TaoViewFromOptions(tao,NULL,"-tao_view"));
187 
188   if (tao->printreason) {
189     if (tao->reason > 0) {
190       PetscCall(PetscPrintf(((PetscObject)tao)->comm,"TAO solve converged due to %s iterations %" PetscInt_FMT "\n",TaoConvergedReasons[tao->reason],tao->niter));
191     } else {
192       PetscCall(PetscPrintf(((PetscObject)tao)->comm,"TAO solve did not converge due to %s iteration %" PetscInt_FMT "\n",TaoConvergedReasons[tao->reason],tao->niter));
193     }
194   }
195   PetscFunctionReturn(0);
196 }
197 
198 /*@
199   TaoSetUp - Sets up the internal data structures for the later use
200   of a Tao solver
201 
202   Collective on tao
203 
204   Input Parameters:
205 . tao - the TAO context
206 
207   Notes:
208   The user will not need to explicitly call TaoSetUp(), as it will
209   automatically be called in TaoSolve().  However, if the user
210   desires to call it explicitly, it should come after TaoCreate()
211   and any TaoSetSomething() routines, but before TaoSolve().
212 
213   Level: advanced
214 
215 .seealso: `TaoCreate()`, `TaoSolve()`
216 @*/
217 PetscErrorCode TaoSetUp(Tao tao)
218 {
219   PetscFunctionBegin;
220   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
221   if (tao->setupcalled) PetscFunctionReturn(0);
222   PetscCall(TaoSetUpEW_Private(tao));
223   PetscCheck(tao->solution,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetSolution");
224   if (tao->ops->setup) PetscCall((*tao->ops->setup)(tao));
225   tao->setupcalled = PETSC_TRUE;
226   PetscFunctionReturn(0);
227 }
228 
229 /*@C
230   TaoDestroy - Destroys the TAO context that was created with
231   TaoCreate()
232 
233   Collective on Tao
234 
235   Input Parameter:
236 . tao - the Tao context
237 
238   Level: beginner
239 
240 .seealso: `TaoCreate()`, `TaoSolve()`
241 @*/
242 PetscErrorCode TaoDestroy(Tao *tao)
243 {
244   PetscFunctionBegin;
245   if (!*tao) PetscFunctionReturn(0);
246   PetscValidHeaderSpecific(*tao,TAO_CLASSID,1);
247   if (--((PetscObject)*tao)->refct > 0) {*tao = NULL;PetscFunctionReturn(0);}
248 
249   if ((*tao)->ops->destroy) {
250     PetscCall((*((*tao))->ops->destroy)(*tao));
251   }
252   PetscCall(KSPDestroy(&(*tao)->ksp));
253   PetscCall(SNESDestroy(&(*tao)->snes_ewdummy));
254   PetscCall(TaoLineSearchDestroy(&(*tao)->linesearch));
255 
256   if ((*tao)->ops->convergencedestroy) {
257     PetscCall((*(*tao)->ops->convergencedestroy)((*tao)->cnvP));
258     if ((*tao)->jacobian_state_inv) {
259       PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
260     }
261   }
262   PetscCall(VecDestroy(&(*tao)->solution));
263   PetscCall(VecDestroy(&(*tao)->gradient));
264   PetscCall(VecDestroy(&(*tao)->ls_res));
265 
266   if ((*tao)->gradient_norm) {
267     PetscCall(PetscObjectDereference((PetscObject)(*tao)->gradient_norm));
268     PetscCall(VecDestroy(&(*tao)->gradient_norm_tmp));
269   }
270 
271   PetscCall(VecDestroy(&(*tao)->XL));
272   PetscCall(VecDestroy(&(*tao)->XU));
273   PetscCall(VecDestroy(&(*tao)->IL));
274   PetscCall(VecDestroy(&(*tao)->IU));
275   PetscCall(VecDestroy(&(*tao)->DE));
276   PetscCall(VecDestroy(&(*tao)->DI));
277   PetscCall(VecDestroy(&(*tao)->constraints));
278   PetscCall(VecDestroy(&(*tao)->constraints_equality));
279   PetscCall(VecDestroy(&(*tao)->constraints_inequality));
280   PetscCall(VecDestroy(&(*tao)->stepdirection));
281   PetscCall(MatDestroy(&(*tao)->hessian_pre));
282   PetscCall(MatDestroy(&(*tao)->hessian));
283   PetscCall(MatDestroy(&(*tao)->ls_jac));
284   PetscCall(MatDestroy(&(*tao)->ls_jac_pre));
285   PetscCall(MatDestroy(&(*tao)->jacobian_pre));
286   PetscCall(MatDestroy(&(*tao)->jacobian));
287   PetscCall(MatDestroy(&(*tao)->jacobian_state_pre));
288   PetscCall(MatDestroy(&(*tao)->jacobian_state));
289   PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
290   PetscCall(MatDestroy(&(*tao)->jacobian_design));
291   PetscCall(MatDestroy(&(*tao)->jacobian_equality));
292   PetscCall(MatDestroy(&(*tao)->jacobian_equality_pre));
293   PetscCall(MatDestroy(&(*tao)->jacobian_inequality));
294   PetscCall(MatDestroy(&(*tao)->jacobian_inequality_pre));
295   PetscCall(ISDestroy(&(*tao)->state_is));
296   PetscCall(ISDestroy(&(*tao)->design_is));
297   PetscCall(VecDestroy(&(*tao)->res_weights_v));
298   PetscCall(TaoCancelMonitors(*tao));
299   if ((*tao)->hist_malloc) {
300     PetscCall(PetscFree4((*tao)->hist_obj,(*tao)->hist_resid,(*tao)->hist_cnorm,(*tao)->hist_lits));
301   }
302   if ((*tao)->res_weights_n) {
303     PetscCall(PetscFree((*tao)->res_weights_rows));
304     PetscCall(PetscFree((*tao)->res_weights_cols));
305     PetscCall(PetscFree((*tao)->res_weights_w));
306   }
307   PetscCall(PetscHeaderDestroy(tao));
308   PetscFunctionReturn(0);
309 }
310 
311 /*@
312    TaoKSPSetUseEW - Sets SNES use Eisenstat-Walker method for
313    computing relative tolerance for linear solvers.
314 
315    Logically Collective on TAO
316 
317    Input Parameters:
318 +  tao - Tao context
319 -  flag - PETSC_TRUE or PETSC_FALSE
320 
321    Notes:
322    See SNESKSPSetUseEW() for customization details.
323 
324    Level: advanced
325 
326    Reference:
327    S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an
328    inexact Newton method", SISC 17 (1), pp.16-32, 1996.
329 
330 .seealso: `SNESKSPSetUseEW()`
331 @*/
332 PetscErrorCode  TaoKSPSetUseEW(Tao tao,PetscBool flag)
333 {
334   PetscFunctionBegin;
335   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
336   PetscValidLogicalCollectiveBool(tao,flag,2);
337   tao->ksp_ewconv = flag;
338   PetscFunctionReturn(0);
339 }
340 
341 /*@
342   TaoSetFromOptions - Sets various Tao parameters from user
343   options.
344 
345   Collective on Tao
346 
347   Input Parameter:
348 . tao - the Tao solver context
349 
350   options Database Keys:
351 + -tao_type <type> - The algorithm that TAO uses (lmvm, nls, etc.)
352 . -tao_gatol <gatol> - absolute error tolerance for ||gradient||
353 . -tao_grtol <grtol> - relative error tolerance for ||gradient||
354 . -tao_gttol <gttol> - reduction of ||gradient|| relative to initial gradient
355 . -tao_max_it <max> - sets maximum number of iterations
356 . -tao_max_funcs <max> - sets maximum number of function evaluations
357 . -tao_fmin <fmin> - stop if function value reaches fmin
358 . -tao_steptol <tol> - stop if trust region radius less than <tol>
359 . -tao_trust0 <t> - initial trust region radius
360 . -tao_monitor - prints function value and residual at each iteration
361 . -tao_smonitor - same as tao_monitor, but truncates very small values
362 . -tao_cmonitor - prints function value, residual, and constraint norm at each iteration
363 . -tao_view_solution - prints solution vector at each iteration
364 . -tao_view_ls_residual - prints least-squares residual vector at each iteration
365 . -tao_view_stepdirection - prints step direction vector at each iteration
366 . -tao_view_gradient - prints gradient vector at each iteration
367 . -tao_draw_solution - graphically view solution vector at each iteration
368 . -tao_draw_step - graphically view step vector at each iteration
369 . -tao_draw_gradient - graphically view gradient at each iteration
370 . -tao_fd_gradient - use gradient computed with finite differences
371 . -tao_fd_hessian - use hessian computed with finite differences
372 . -tao_mf_hessian - use matrix-free hessian computed with finite differences
373 . -tao_cancelmonitors - cancels all monitors (except those set with command line)
374 . -tao_view - prints information about the Tao after solving
375 - -tao_converged_reason - prints the reason TAO stopped iterating
376 
377   Notes:
378   To see all options, run your program with the -help option or consult the
379   user's manual. Should be called after TaoCreate() but before TaoSolve()
380 
381   Level: beginner
382 @*/
383 PetscErrorCode TaoSetFromOptions(Tao tao)
384 {
385   TaoType        default_type = TAOLMVM;
386   char           type[256], monfilename[PETSC_MAX_PATH_LEN];
387   PetscViewer    monviewer;
388   PetscBool      flg;
389   MPI_Comm       comm;
390 
391   PetscFunctionBegin;
392   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
393   PetscCall(PetscObjectGetComm((PetscObject)tao,&comm));
394 
395   /* So no warnings are given about unused options */
396   PetscCall(PetscOptionsHasName(((PetscObject)tao)->options,((PetscObject)tao)->prefix,"-tao_ls_type",&flg));
397 
398   PetscObjectOptionsBegin((PetscObject)tao);
399   {
400     if (((PetscObject)tao)->type_name) default_type = ((PetscObject)tao)->type_name;
401     /* Check for type from options */
402     PetscCall(PetscOptionsFList("-tao_type","Tao Solver type","TaoSetType",TaoList,default_type,type,256,&flg));
403     if (flg) {
404       PetscCall(TaoSetType(tao,type));
405     } else if (!((PetscObject)tao)->type_name) {
406       PetscCall(TaoSetType(tao,default_type));
407     }
408 
409     PetscCall(PetscOptionsReal("-tao_catol","Stop if constraints violations within","TaoSetConstraintTolerances",tao->catol,&tao->catol,&flg));
410     if (flg) tao->catol_changed = PETSC_TRUE;
411     PetscCall(PetscOptionsReal("-tao_crtol","Stop if relative constraint violations within","TaoSetConstraintTolerances",tao->crtol,&tao->crtol,&flg));
412     if (flg) tao->crtol_changed = PETSC_TRUE;
413     PetscCall(PetscOptionsReal("-tao_gatol","Stop if norm of gradient less than","TaoSetTolerances",tao->gatol,&tao->gatol,&flg));
414     if (flg) tao->gatol_changed = PETSC_TRUE;
415     PetscCall(PetscOptionsReal("-tao_grtol","Stop if norm of gradient divided by the function value is less than","TaoSetTolerances",tao->grtol,&tao->grtol,&flg));
416     if (flg) tao->grtol_changed = PETSC_TRUE;
417     PetscCall(PetscOptionsReal("-tao_gttol","Stop if the norm of the gradient is less than the norm of the initial gradient times tol","TaoSetTolerances",tao->gttol,&tao->gttol,&flg));
418     if (flg) tao->gttol_changed = PETSC_TRUE;
419     PetscCall(PetscOptionsInt("-tao_max_it","Stop if iteration number exceeds","TaoSetMaximumIterations",tao->max_it,&tao->max_it,&flg));
420     if (flg) tao->max_it_changed = PETSC_TRUE;
421     PetscCall(PetscOptionsInt("-tao_max_funcs","Stop if number of function evaluations exceeds","TaoSetMaximumFunctionEvaluations",tao->max_funcs,&tao->max_funcs,&flg));
422     if (flg) tao->max_funcs_changed = PETSC_TRUE;
423     PetscCall(PetscOptionsReal("-tao_fmin","Stop if function less than","TaoSetFunctionLowerBound",tao->fmin,&tao->fmin,&flg));
424     if (flg) tao->fmin_changed = PETSC_TRUE;
425     PetscCall(PetscOptionsReal("-tao_steptol","Stop if step size or trust region radius less than","",tao->steptol,&tao->steptol,&flg));
426     if (flg) tao->steptol_changed = PETSC_TRUE;
427     PetscCall(PetscOptionsReal("-tao_trust0","Initial trust region radius","TaoSetTrustRegionRadius",tao->trust0,&tao->trust0,&flg));
428     if (flg) tao->trust0_changed = PETSC_TRUE;
429     PetscCall(PetscOptionsString("-tao_view_solution","view solution vector after each evaluation","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
430     if (flg) {
431       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
432       PetscCall(TaoSetMonitor(tao,TaoSolutionMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
433     }
434 
435     PetscCall(PetscOptionsBool("-tao_converged_reason","Print reason for TAO converged","TaoSolve",tao->printreason,&tao->printreason,NULL));
436     PetscCall(PetscOptionsString("-tao_view_gradient","view gradient vector after each evaluation","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
437     if (flg) {
438       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
439       PetscCall(TaoSetMonitor(tao,TaoGradientMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
440     }
441 
442     PetscCall(PetscOptionsString("-tao_view_stepdirection","view step direction vector after each iteration","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
443     if (flg) {
444       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
445       PetscCall(TaoSetMonitor(tao,TaoStepDirectionMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
446     }
447 
448     PetscCall(PetscOptionsString("-tao_view_residual","view least-squares residual vector after each evaluation","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
449     if (flg) {
450       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
451       PetscCall(TaoSetMonitor(tao,TaoResidualMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
452     }
453 
454     PetscCall(PetscOptionsString("-tao_monitor","Use the default convergence monitor","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
455     if (flg) {
456       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
457       PetscCall(TaoSetMonitor(tao,TaoMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
458     }
459 
460     PetscCall(PetscOptionsString("-tao_gmonitor","Use the convergence monitor with extra globalization info","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
461     if (flg) {
462       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
463       PetscCall(TaoSetMonitor(tao,TaoDefaultGMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
464     }
465 
466     PetscCall(PetscOptionsString("-tao_smonitor","Use the short convergence monitor","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
467     if (flg) {
468       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
469       PetscCall(TaoSetMonitor(tao,TaoDefaultSMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
470     }
471 
472     PetscCall(PetscOptionsString("-tao_cmonitor","Use the default convergence monitor with constraint norm","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
473     if (flg) {
474       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
475       PetscCall(TaoSetMonitor(tao,TaoDefaultCMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
476     }
477 
478     flg = PETSC_FALSE;
479     PetscCall(PetscOptionsBool("-tao_cancelmonitors","cancel all monitors and call any registered destroy routines","TaoCancelMonitors",flg,&flg,NULL));
480     if (flg) PetscCall(TaoCancelMonitors(tao));
481 
482     flg = PETSC_FALSE;
483     PetscCall(PetscOptionsBool("-tao_draw_solution","Plot solution vector at each iteration","TaoSetMonitor",flg,&flg,NULL));
484     if (flg) {
485       TaoMonitorDrawCtx drawctx;
486       PetscInt          howoften = 1;
487       PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&drawctx));
488       PetscCall(TaoSetMonitor(tao,TaoDrawSolutionMonitor,drawctx,(PetscErrorCode (*)(void**))TaoMonitorDrawCtxDestroy));
489     }
490 
491     flg = PETSC_FALSE;
492     PetscCall(PetscOptionsBool("-tao_draw_step","plots step direction at each iteration","TaoSetMonitor",flg,&flg,NULL));
493     if (flg) {
494       PetscCall(TaoSetMonitor(tao,TaoDrawStepMonitor,NULL,NULL));
495     }
496 
497     flg = PETSC_FALSE;
498     PetscCall(PetscOptionsBool("-tao_draw_gradient","plots gradient at each iteration","TaoSetMonitor",flg,&flg,NULL));
499     if (flg) {
500       TaoMonitorDrawCtx drawctx;
501       PetscInt          howoften = 1;
502       PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&drawctx));
503       PetscCall(TaoSetMonitor(tao,TaoDrawGradientMonitor,drawctx,(PetscErrorCode (*)(void**))TaoMonitorDrawCtxDestroy));
504     }
505     flg = PETSC_FALSE;
506     PetscCall(PetscOptionsBool("-tao_fd_gradient","compute gradient using finite differences","TaoDefaultComputeGradient",flg,&flg,NULL));
507     if (flg) {
508       PetscCall(TaoSetGradient(tao,NULL,TaoDefaultComputeGradient,NULL));
509     }
510     flg = PETSC_FALSE;
511     PetscCall(PetscOptionsBool("-tao_fd_hessian","compute hessian using finite differences","TaoDefaultComputeHessian",flg,&flg,NULL));
512     if (flg) {
513       Mat H;
514 
515       PetscCall(MatCreate(PetscObjectComm((PetscObject)tao),&H));
516       PetscCall(MatSetType(H,MATAIJ));
517       PetscCall(TaoSetHessian(tao,H,H,TaoDefaultComputeHessian,NULL));
518       PetscCall(MatDestroy(&H));
519     }
520     flg = PETSC_FALSE;
521     PetscCall(PetscOptionsBool("-tao_mf_hessian","compute matrix-free hessian using finite differences","TaoDefaultComputeHessianMFFD",flg,&flg,NULL));
522     if (flg) {
523       Mat H;
524 
525       PetscCall(MatCreate(PetscObjectComm((PetscObject)tao),&H));
526       PetscCall(TaoSetHessian(tao,H,H,TaoDefaultComputeHessianMFFD,NULL));
527       PetscCall(MatDestroy(&H));
528     }
529     flg = PETSC_FALSE;
530     PetscCall(PetscOptionsBool("-tao_recycle_history","enable recycling/re-using information from the previous TaoSolve() call for some algorithms","TaoSetRecycleHistory",flg,&flg,NULL));
531     if (flg) {
532       PetscCall(TaoSetRecycleHistory(tao,PETSC_TRUE));
533     }
534     PetscCall(PetscOptionsEnum("-tao_subset_type","subset type","",TaoSubSetTypes,(PetscEnum)tao->subset_type,(PetscEnum*)&tao->subset_type,NULL));
535 
536     if (tao->ksp) {
537       PetscCall(PetscOptionsBool("-tao_ksp_ew","Use Eisentat-Walker linear system convergence test","TaoKSPSetUseEW",tao->ksp_ewconv,&tao->ksp_ewconv,NULL));
538       PetscCall(TaoKSPSetUseEW(tao,tao->ksp_ewconv));
539     }
540 
541     if (tao->linesearch) {
542       PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
543     }
544 
545     if (tao->ops->setfromoptions) {
546       PetscCall((*tao->ops->setfromoptions)(PetscOptionsObject,tao));
547     }
548   }
549   PetscOptionsEnd();
550   PetscFunctionReturn(0);
551 }
552 
553 /*@C
554    TaoViewFromOptions - View from Options
555 
556    Collective on Tao
557 
558    Input Parameters:
559 +  A - the  Tao context
560 .  obj - Optional object
561 -  name - command line option
562 
563    Level: intermediate
564 .seealso: `Tao`, `TaoView`, `PetscObjectViewFromOptions()`, `TaoCreate()`
565 @*/
566 PetscErrorCode  TaoViewFromOptions(Tao A,PetscObject obj,const char name[])
567 {
568   PetscFunctionBegin;
569   PetscValidHeaderSpecific(A,TAO_CLASSID,1);
570   PetscCall(PetscObjectViewFromOptions((PetscObject)A,obj,name));
571   PetscFunctionReturn(0);
572 }
573 
574 /*@C
575   TaoView - Prints information about the Tao
576 
577   Collective on Tao
578 
579   InputParameters:
580 + tao - the Tao context
581 - viewer - visualization context
582 
583   Options Database Key:
584 . -tao_view - Calls TaoView() at the end of TaoSolve()
585 
586   Notes:
587   The available visualization contexts include
588 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
589 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
590          output where only the first processor opens
591          the file.  All other processors send their
592          data to the first processor to print.
593 
594   Level: beginner
595 
596 .seealso: `PetscViewerASCIIOpen()`
597 @*/
598 PetscErrorCode TaoView(Tao tao, PetscViewer viewer)
599 {
600   PetscBool           isascii,isstring;
601   TaoType             type;
602 
603   PetscFunctionBegin;
604   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
605   if (!viewer) {
606     PetscCall(PetscViewerASCIIGetStdout(((PetscObject)tao)->comm,&viewer));
607   }
608   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
609   PetscCheckSameComm(tao,1,viewer,2);
610 
611   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
612   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring));
613   if (isascii) {
614     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)tao,viewer));
615 
616     if (tao->ops->view) {
617       PetscCall(PetscViewerASCIIPushTab(viewer));
618       PetscCall((*tao->ops->view)(tao,viewer));
619       PetscCall(PetscViewerASCIIPopTab(viewer));
620     }
621     if (tao->linesearch) {
622       PetscCall(PetscViewerASCIIPushTab(viewer));
623       PetscCall(TaoLineSearchView(tao->linesearch,viewer));
624       PetscCall(PetscViewerASCIIPopTab(viewer));
625     }
626     if (tao->ksp) {
627       PetscCall(PetscViewerASCIIPushTab(viewer));
628       PetscCall(KSPView(tao->ksp,viewer));
629       PetscCall(PetscViewerASCIIPrintf(viewer,"total KSP iterations: %" PetscInt_FMT "\n",tao->ksp_tot_its));
630       PetscCall(PetscViewerASCIIPopTab(viewer));
631     }
632 
633     PetscCall(PetscViewerASCIIPushTab(viewer));
634 
635     if (tao->XL || tao->XU) {
636       PetscCall(PetscViewerASCIIPrintf(viewer,"Active Set subset type: %s\n",TaoSubSetTypes[tao->subset_type]));
637     }
638 
639     PetscCall(PetscViewerASCIIPrintf(viewer,"convergence tolerances: gatol=%g,",(double)tao->gatol));
640     PetscCall(PetscViewerASCIIPrintf(viewer," steptol=%g,",(double)tao->steptol));
641     PetscCall(PetscViewerASCIIPrintf(viewer," gttol=%g\n",(double)tao->gttol));
642     PetscCall(PetscViewerASCIIPrintf(viewer,"Residual in Function/Gradient:=%g\n",(double)tao->residual));
643 
644     if (tao->constrained) {
645       PetscCall(PetscViewerASCIIPrintf(viewer,"convergence tolerances:"));
646       PetscCall(PetscViewerASCIIPrintf(viewer," catol=%g,",(double)tao->catol));
647       PetscCall(PetscViewerASCIIPrintf(viewer," crtol=%g\n",(double)tao->crtol));
648       PetscCall(PetscViewerASCIIPrintf(viewer,"Residual in Constraints:=%g\n",(double)tao->cnorm));
649     }
650 
651     if (tao->trust < tao->steptol) {
652       PetscCall(PetscViewerASCIIPrintf(viewer,"convergence tolerances: steptol=%g\n",(double)tao->steptol));
653       PetscCall(PetscViewerASCIIPrintf(viewer,"Final trust region radius:=%g\n",(double)tao->trust));
654     }
655 
656     if (tao->fmin>-1.e25) {
657       PetscCall(PetscViewerASCIIPrintf(viewer,"convergence tolerances: function minimum=%g\n",(double)tao->fmin));
658     }
659     PetscCall(PetscViewerASCIIPrintf(viewer,"Objective value=%g\n",(double)tao->fc));
660 
661     PetscCall(PetscViewerASCIIPrintf(viewer,"total number of iterations=%" PetscInt_FMT ",          ",tao->niter));
662     PetscCall(PetscViewerASCIIPrintf(viewer,"              (max: %" PetscInt_FMT ")\n",tao->max_it));
663 
664     if (tao->nfuncs>0) {
665       PetscCall(PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%" PetscInt_FMT ",",tao->nfuncs));
666       PetscCall(PetscViewerASCIIPrintf(viewer,"                max: %" PetscInt_FMT "\n",tao->max_funcs));
667     }
668     if (tao->ngrads>0) {
669       PetscCall(PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%" PetscInt_FMT ",",tao->ngrads));
670       PetscCall(PetscViewerASCIIPrintf(viewer,"                max: %" PetscInt_FMT "\n",tao->max_funcs));
671     }
672     if (tao->nfuncgrads>0) {
673       PetscCall(PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%" PetscInt_FMT ",",tao->nfuncgrads));
674       PetscCall(PetscViewerASCIIPrintf(viewer,"    (max: %" PetscInt_FMT ")\n",tao->max_funcs));
675     }
676     if (tao->nhess>0) {
677       PetscCall(PetscViewerASCIIPrintf(viewer,"total number of Hessian evaluations=%" PetscInt_FMT "\n",tao->nhess));
678     }
679     if (tao->nconstraints>0) {
680       PetscCall(PetscViewerASCIIPrintf(viewer,"total number of constraint function evaluations=%" PetscInt_FMT "\n",tao->nconstraints));
681     }
682     if (tao->njac>0) {
683       PetscCall(PetscViewerASCIIPrintf(viewer,"total number of Jacobian evaluations=%" PetscInt_FMT "\n",tao->njac));
684     }
685 
686     if (tao->reason>0) {
687       PetscCall(PetscViewerASCIIPrintf(viewer,    "Solution converged: "));
688       switch (tao->reason) {
689       case TAO_CONVERGED_GATOL:
690         PetscCall(PetscViewerASCIIPrintf(viewer," ||g(X)|| <= gatol\n"));
691         break;
692       case TAO_CONVERGED_GRTOL:
693         PetscCall(PetscViewerASCIIPrintf(viewer," ||g(X)||/|f(X)| <= grtol\n"));
694         break;
695       case TAO_CONVERGED_GTTOL:
696         PetscCall(PetscViewerASCIIPrintf(viewer," ||g(X)||/||g(X0)|| <= gttol\n"));
697         break;
698       case TAO_CONVERGED_STEPTOL:
699         PetscCall(PetscViewerASCIIPrintf(viewer," Steptol -- step size small\n"));
700         break;
701       case TAO_CONVERGED_MINF:
702         PetscCall(PetscViewerASCIIPrintf(viewer," Minf --  f < fmin\n"));
703         break;
704       case TAO_CONVERGED_USER:
705         PetscCall(PetscViewerASCIIPrintf(viewer," User Terminated\n"));
706         break;
707       default:
708         PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
709         break;
710       }
711     } else {
712       PetscCall(PetscViewerASCIIPrintf(viewer,"Solver terminated: %d",tao->reason));
713       switch (tao->reason) {
714       case TAO_DIVERGED_MAXITS:
715         PetscCall(PetscViewerASCIIPrintf(viewer," Maximum Iterations\n"));
716         break;
717       case TAO_DIVERGED_NAN:
718         PetscCall(PetscViewerASCIIPrintf(viewer," NAN or Inf encountered\n"));
719         break;
720       case TAO_DIVERGED_MAXFCN:
721         PetscCall(PetscViewerASCIIPrintf(viewer," Maximum Function Evaluations\n"));
722         break;
723       case TAO_DIVERGED_LS_FAILURE:
724         PetscCall(PetscViewerASCIIPrintf(viewer," Line Search Failure\n"));
725         break;
726       case TAO_DIVERGED_TR_REDUCTION:
727         PetscCall(PetscViewerASCIIPrintf(viewer," Trust Region too small\n"));
728         break;
729       case TAO_DIVERGED_USER:
730         PetscCall(PetscViewerASCIIPrintf(viewer," User Terminated\n"));
731         break;
732       default:
733         PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
734         break;
735       }
736     }
737     PetscCall(PetscViewerASCIIPopTab(viewer));
738   } else if (isstring) {
739     PetscCall(TaoGetType(tao,&type));
740     PetscCall(PetscViewerStringSPrintf(viewer," %-3.3s",type));
741   }
742   PetscFunctionReturn(0);
743 }
744 
745 /*@
746   TaoSetRecycleHistory - Sets the boolean flag to enable/disable re-using
747   iterate information from the previous TaoSolve(). This feature is disabled by
748   default.
749 
750   For conjugate gradient methods (BNCG), this re-uses the latest search direction
751   from the previous TaoSolve() call when computing the first search direction in a
752   new solution. By default, CG methods set the first search direction to the
753   negative gradient.
754 
755   For quasi-Newton family of methods (BQNLS, BQNKLS, BQNKTR, BQNKTL), this re-uses
756   the accumulated quasi-Newton Hessian approximation from the previous TaoSolve()
757   call. By default, QN family of methods reset the initial Hessian approximation to
758   the identity matrix.
759 
760   For any other algorithm, this setting has no effect.
761 
762   Logically collective on Tao
763 
764   Input Parameters:
765 + tao - the Tao context
766 - recycle - boolean flag
767 
768   Options Database Keys:
769 . -tao_recycle_history <true,false> - reuse the history
770 
771   Level: intermediate
772 
773 .seealso: `TaoSetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
774 
775 @*/
776 PetscErrorCode TaoSetRecycleHistory(Tao tao, PetscBool recycle)
777 {
778   PetscFunctionBegin;
779   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
780   PetscValidLogicalCollectiveBool(tao,recycle,2);
781   tao->recycle = recycle;
782   PetscFunctionReturn(0);
783 }
784 
785 /*@
786   TaoGetRecycleHistory - Retrieve the boolean flag for re-using iterate information
787   from the previous TaoSolve(). This feature is disabled by default.
788 
789   Logically collective on Tao
790 
791   Input Parameters:
792 . tao - the Tao context
793 
794   Output Parameters:
795 . recycle - boolean flag
796 
797   Level: intermediate
798 
799 .seealso: `TaoGetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
800 
801 @*/
802 PetscErrorCode TaoGetRecycleHistory(Tao tao, PetscBool *recycle)
803 {
804   PetscFunctionBegin;
805   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
806   PetscValidBoolPointer(recycle,2);
807   *recycle = tao->recycle;
808   PetscFunctionReturn(0);
809 }
810 
811 /*@
812   TaoSetTolerances - Sets parameters used in TAO convergence tests
813 
814   Logically collective on Tao
815 
816   Input Parameters:
817 + tao - the Tao context
818 . gatol - stop if norm of gradient is less than this
819 . grtol - stop if relative norm of gradient is less than this
820 - gttol - stop if norm of gradient is reduced by this factor
821 
822   Options Database Keys:
823 + -tao_gatol <gatol> - Sets gatol
824 . -tao_grtol <grtol> - Sets grtol
825 - -tao_gttol <gttol> - Sets gttol
826 
827   Stopping Criteria:
828 $ ||g(X)||                            <= gatol
829 $ ||g(X)|| / |f(X)|                   <= grtol
830 $ ||g(X)|| / ||g(X0)||                <= gttol
831 
832   Notes:
833   Use PETSC_DEFAULT to leave one or more tolerances unchanged.
834 
835   Level: beginner
836 
837 .seealso: `TaoGetTolerances()`
838 
839 @*/
840 PetscErrorCode TaoSetTolerances(Tao tao, PetscReal gatol, PetscReal grtol, PetscReal gttol)
841 {
842   PetscFunctionBegin;
843   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
844   PetscValidLogicalCollectiveReal(tao,gatol,2);
845   PetscValidLogicalCollectiveReal(tao,grtol,3);
846   PetscValidLogicalCollectiveReal(tao,gttol,4);
847 
848   if (gatol != PETSC_DEFAULT) {
849     if (gatol<0) {
850       PetscCall(PetscInfo(tao,"Tried to set negative gatol -- ignored.\n"));
851     } else {
852       tao->gatol = PetscMax(0,gatol);
853       tao->gatol_changed = PETSC_TRUE;
854     }
855   }
856 
857   if (grtol != PETSC_DEFAULT) {
858     if (grtol<0) {
859       PetscCall(PetscInfo(tao,"Tried to set negative grtol -- ignored.\n"));
860     } else {
861       tao->grtol = PetscMax(0,grtol);
862       tao->grtol_changed = PETSC_TRUE;
863     }
864   }
865 
866   if (gttol != PETSC_DEFAULT) {
867     if (gttol<0) {
868       PetscCall(PetscInfo(tao,"Tried to set negative gttol -- ignored.\n"));
869     } else {
870       tao->gttol = PetscMax(0,gttol);
871       tao->gttol_changed = PETSC_TRUE;
872     }
873   }
874   PetscFunctionReturn(0);
875 }
876 
877 /*@
878   TaoSetConstraintTolerances - Sets constraint tolerance parameters used in TAO  convergence tests
879 
880   Logically collective on Tao
881 
882   Input Parameters:
883 + tao - the Tao context
884 . catol - absolute constraint tolerance, constraint norm must be less than catol for used for gatol convergence criteria
885 - crtol - relative constraint tolerance, constraint norm must be less than crtol for used for gatol, gttol convergence criteria
886 
887   Options Database Keys:
888 + -tao_catol <catol> - Sets catol
889 - -tao_crtol <crtol> - Sets crtol
890 
891   Notes:
892   Use PETSC_DEFAULT to leave any tolerance unchanged.
893 
894   Level: intermediate
895 
896 .seealso: `TaoGetTolerances()`, `TaoGetConstraintTolerances()`, `TaoSetTolerances()`
897 
898 @*/
899 PetscErrorCode TaoSetConstraintTolerances(Tao tao, PetscReal catol, PetscReal crtol)
900 {
901   PetscFunctionBegin;
902   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
903   PetscValidLogicalCollectiveReal(tao,catol,2);
904   PetscValidLogicalCollectiveReal(tao,crtol,3);
905 
906   if (catol != PETSC_DEFAULT) {
907     if (catol<0) {
908       PetscCall(PetscInfo(tao,"Tried to set negative catol -- ignored.\n"));
909     } else {
910       tao->catol = PetscMax(0,catol);
911       tao->catol_changed = PETSC_TRUE;
912     }
913   }
914 
915   if (crtol != PETSC_DEFAULT) {
916     if (crtol<0) {
917       PetscCall(PetscInfo(tao,"Tried to set negative crtol -- ignored.\n"));
918     } else {
919       tao->crtol = PetscMax(0,crtol);
920       tao->crtol_changed = PETSC_TRUE;
921     }
922   }
923   PetscFunctionReturn(0);
924 }
925 
926 /*@
927   TaoGetConstraintTolerances - Gets constraint tolerance parameters used in TAO  convergence tests
928 
929   Not ollective
930 
931   Input Parameter:
932 . tao - the Tao context
933 
934   Output Parameters:
935 + catol - absolute constraint tolerance, constraint norm must be less than catol for used for gatol convergence criteria
936 - crtol - relative constraint tolerance, constraint norm must be less than crtol for used for gatol, gttol convergence criteria
937 
938   Level: intermediate
939 
940 .seealso: `TaoGetTolerances()`, `TaoSetTolerances()`, `TaoSetConstraintTolerances()`
941 
942 @*/
943 PetscErrorCode TaoGetConstraintTolerances(Tao tao, PetscReal *catol, PetscReal *crtol)
944 {
945   PetscFunctionBegin;
946   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
947   if (catol) *catol = tao->catol;
948   if (crtol) *crtol = tao->crtol;
949   PetscFunctionReturn(0);
950 }
951 
952 /*@
953    TaoSetFunctionLowerBound - Sets a bound on the solution objective value.
954    When an approximate solution with an objective value below this number
955    has been found, the solver will terminate.
956 
957    Logically Collective on Tao
958 
959    Input Parameters:
960 +  tao - the Tao solver context
961 -  fmin - the tolerance
962 
963    Options Database Keys:
964 .    -tao_fmin <fmin> - sets the minimum function value
965 
966    Level: intermediate
967 
968 .seealso: `TaoSetTolerances()`
969 @*/
970 PetscErrorCode TaoSetFunctionLowerBound(Tao tao,PetscReal fmin)
971 {
972   PetscFunctionBegin;
973   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
974   PetscValidLogicalCollectiveReal(tao,fmin,2);
975   tao->fmin = fmin;
976   tao->fmin_changed = PETSC_TRUE;
977   PetscFunctionReturn(0);
978 }
979 
980 /*@
981    TaoGetFunctionLowerBound - Gets the bound on the solution objective value.
982    When an approximate solution with an objective value below this number
983    has been found, the solver will terminate.
984 
985    Not collective on Tao
986 
987    Input Parameters:
988 .  tao - the Tao solver context
989 
990    OutputParameters:
991 .  fmin - the minimum function value
992 
993    Level: intermediate
994 
995 .seealso: `TaoSetFunctionLowerBound()`
996 @*/
997 PetscErrorCode TaoGetFunctionLowerBound(Tao tao,PetscReal *fmin)
998 {
999   PetscFunctionBegin;
1000   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1001   PetscValidRealPointer(fmin,2);
1002   *fmin = tao->fmin;
1003   PetscFunctionReturn(0);
1004 }
1005 
1006 /*@
1007    TaoSetMaximumFunctionEvaluations - Sets a maximum number of
1008    function evaluations.
1009 
1010    Logically Collective on Tao
1011 
1012    Input Parameters:
1013 +  tao - the Tao solver context
1014 -  nfcn - the maximum number of function evaluations (>=0)
1015 
1016    Options Database Keys:
1017 .    -tao_max_funcs <nfcn> - sets the maximum number of function evaluations
1018 
1019    Level: intermediate
1020 
1021 .seealso: `TaoSetTolerances()`, `TaoSetMaximumIterations()`
1022 @*/
1023 
1024 PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao tao,PetscInt nfcn)
1025 {
1026   PetscFunctionBegin;
1027   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1028   PetscValidLogicalCollectiveInt(tao,nfcn,2);
1029   if (nfcn >= 0) { tao->max_funcs = PetscMax(0,nfcn); }
1030   else { tao->max_funcs = -1; }
1031   tao->max_funcs_changed = PETSC_TRUE;
1032   PetscFunctionReturn(0);
1033 }
1034 
1035 /*@
1036    TaoGetMaximumFunctionEvaluations - Sets a maximum number of
1037    function evaluations.
1038 
1039    Not Collective
1040 
1041    Input Parameters:
1042 .  tao - the Tao solver context
1043 
1044    Output Parameters:
1045 .  nfcn - the maximum number of function evaluations
1046 
1047    Level: intermediate
1048 
1049 .seealso: `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1050 @*/
1051 
1052 PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao tao,PetscInt *nfcn)
1053 {
1054   PetscFunctionBegin;
1055   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1056   PetscValidIntPointer(nfcn,2);
1057   *nfcn = tao->max_funcs;
1058   PetscFunctionReturn(0);
1059 }
1060 
1061 /*@
1062    TaoGetCurrentFunctionEvaluations - Get current number of
1063    function evaluations.
1064 
1065    Not Collective
1066 
1067    Input Parameters:
1068 .  tao - the Tao solver context
1069 
1070    Output Parameters:
1071 .  nfuncs - the current number of function evaluations (maximum between gradient and function evaluations)
1072 
1073    Level: intermediate
1074 
1075 .seealso: `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1076 @*/
1077 
1078 PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao tao,PetscInt *nfuncs)
1079 {
1080   PetscFunctionBegin;
1081   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1082   PetscValidIntPointer(nfuncs,2);
1083   *nfuncs = PetscMax(tao->nfuncs,tao->nfuncgrads);
1084   PetscFunctionReturn(0);
1085 }
1086 
1087 /*@
1088    TaoSetMaximumIterations - Sets a maximum number of iterates.
1089 
1090    Logically Collective on Tao
1091 
1092    Input Parameters:
1093 +  tao - the Tao solver context
1094 -  maxits - the maximum number of iterates (>=0)
1095 
1096    Options Database Keys:
1097 .    -tao_max_it <its> - sets the maximum number of iterations
1098 
1099    Level: intermediate
1100 
1101 .seealso: `TaoSetTolerances()`, `TaoSetMaximumFunctionEvaluations()`
1102 @*/
1103 PetscErrorCode TaoSetMaximumIterations(Tao tao,PetscInt maxits)
1104 {
1105   PetscFunctionBegin;
1106   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1107   PetscValidLogicalCollectiveInt(tao,maxits,2);
1108   tao->max_it = PetscMax(0,maxits);
1109   tao->max_it_changed = PETSC_TRUE;
1110   PetscFunctionReturn(0);
1111 }
1112 
1113 /*@
1114    TaoGetMaximumIterations - Sets a maximum number of iterates.
1115 
1116    Not Collective
1117 
1118    Input Parameters:
1119 .  tao - the Tao solver context
1120 
1121    Output Parameters:
1122 .  maxits - the maximum number of iterates
1123 
1124    Level: intermediate
1125 
1126 .seealso: `TaoSetMaximumIterations()`, `TaoGetMaximumFunctionEvaluations()`
1127 @*/
1128 PetscErrorCode TaoGetMaximumIterations(Tao tao,PetscInt *maxits)
1129 {
1130   PetscFunctionBegin;
1131   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1132   PetscValidIntPointer(maxits,2);
1133   *maxits = tao->max_it;
1134   PetscFunctionReturn(0);
1135 }
1136 
1137 /*@
1138    TaoSetInitialTrustRegionRadius - Sets the initial trust region radius.
1139 
1140    Logically collective on Tao
1141 
1142    Input Parameters:
1143 +  tao - a TAO optimization solver
1144 -  radius - the trust region radius
1145 
1146    Level: intermediate
1147 
1148    Options Database Key:
1149 .  -tao_trust0 <t0> - sets initial trust region radius
1150 
1151 .seealso: `TaoGetTrustRegionRadius()`, `TaoSetTrustRegionTolerance()`
1152 @*/
1153 PetscErrorCode TaoSetInitialTrustRegionRadius(Tao tao, PetscReal radius)
1154 {
1155   PetscFunctionBegin;
1156   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1157   PetscValidLogicalCollectiveReal(tao,radius,2);
1158   tao->trust0 = PetscMax(0.0,radius);
1159   tao->trust0_changed = PETSC_TRUE;
1160   PetscFunctionReturn(0);
1161 }
1162 
1163 /*@
1164    TaoGetInitialTrustRegionRadius - Sets the initial trust region radius.
1165 
1166    Not Collective
1167 
1168    Input Parameter:
1169 .  tao - a TAO optimization solver
1170 
1171    Output Parameter:
1172 .  radius - the trust region radius
1173 
1174    Level: intermediate
1175 
1176 .seealso: `TaoSetInitialTrustRegionRadius()`, `TaoGetCurrentTrustRegionRadius()`
1177 @*/
1178 PetscErrorCode TaoGetInitialTrustRegionRadius(Tao tao, PetscReal *radius)
1179 {
1180   PetscFunctionBegin;
1181   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1182   PetscValidRealPointer(radius,2);
1183   *radius = tao->trust0;
1184   PetscFunctionReturn(0);
1185 }
1186 
1187 /*@
1188    TaoGetCurrentTrustRegionRadius - Gets the current trust region radius.
1189 
1190    Not Collective
1191 
1192    Input Parameter:
1193 .  tao - a TAO optimization solver
1194 
1195    Output Parameter:
1196 .  radius - the trust region radius
1197 
1198    Level: intermediate
1199 
1200 .seealso: `TaoSetInitialTrustRegionRadius()`, `TaoGetInitialTrustRegionRadius()`
1201 @*/
1202 PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao tao, PetscReal *radius)
1203 {
1204   PetscFunctionBegin;
1205   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1206   PetscValidRealPointer(radius,2);
1207   *radius = tao->trust;
1208   PetscFunctionReturn(0);
1209 }
1210 
1211 /*@
1212   TaoGetTolerances - gets the current values of tolerances
1213 
1214   Not Collective
1215 
1216   Input Parameter:
1217 . tao - the Tao context
1218 
1219   Output Parameters:
1220 + gatol - stop if norm of gradient is less than this
1221 . grtol - stop if relative norm of gradient is less than this
1222 - gttol - stop if norm of gradient is reduced by a this factor
1223 
1224   Note: NULL can be used as an argument if not all tolerances values are needed
1225 
1226 .seealso `TaoSetTolerances()`
1227 
1228   Level: intermediate
1229 @*/
1230 PetscErrorCode TaoGetTolerances(Tao tao, PetscReal *gatol, PetscReal *grtol, PetscReal *gttol)
1231 {
1232   PetscFunctionBegin;
1233   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1234   if (gatol) *gatol = tao->gatol;
1235   if (grtol) *grtol = tao->grtol;
1236   if (gttol) *gttol = tao->gttol;
1237   PetscFunctionReturn(0);
1238 }
1239 
1240 /*@
1241   TaoGetKSP - Gets the linear solver used by the optimization solver.
1242   Application writers should use TaoGetKSP if they need direct access
1243   to the PETSc KSP object.
1244 
1245   Not Collective
1246 
1247    Input Parameters:
1248 .  tao - the TAO solver
1249 
1250    Output Parameters:
1251 .  ksp - the KSP linear solver used in the optimization solver
1252 
1253    Level: intermediate
1254 
1255 @*/
1256 PetscErrorCode TaoGetKSP(Tao tao, KSP *ksp)
1257 {
1258   PetscFunctionBegin;
1259   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1260   PetscValidPointer(ksp,2);
1261   *ksp = tao->ksp;
1262   PetscFunctionReturn(0);
1263 }
1264 
1265 /*@
1266    TaoGetLinearSolveIterations - Gets the total number of linear iterations
1267    used by the TAO solver
1268 
1269    Not Collective
1270 
1271    Input Parameter:
1272 .  tao - TAO context
1273 
1274    Output Parameter:
1275 .  lits - number of linear iterations
1276 
1277    Notes:
1278    This counter is reset to zero for each successive call to TaoSolve()
1279 
1280    Level: intermediate
1281 
1282 .seealso: `TaoGetKSP()`
1283 @*/
1284 PetscErrorCode TaoGetLinearSolveIterations(Tao tao, PetscInt *lits)
1285 {
1286   PetscFunctionBegin;
1287   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1288   PetscValidIntPointer(lits,2);
1289   *lits = tao->ksp_tot_its;
1290   PetscFunctionReturn(0);
1291 }
1292 
1293 /*@
1294   TaoGetLineSearch - Gets the line search used by the optimization solver.
1295   Application writers should use TaoGetLineSearch if they need direct access
1296   to the TaoLineSearch object.
1297 
1298   Not Collective
1299 
1300    Input Parameters:
1301 .  tao - the TAO solver
1302 
1303    Output Parameters:
1304 .  ls - the line search used in the optimization solver
1305 
1306    Level: intermediate
1307 
1308 @*/
1309 PetscErrorCode TaoGetLineSearch(Tao tao, TaoLineSearch *ls)
1310 {
1311   PetscFunctionBegin;
1312   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1313   PetscValidPointer(ls,2);
1314   *ls = tao->linesearch;
1315   PetscFunctionReturn(0);
1316 }
1317 
1318 /*@
1319   TaoAddLineSearchCounts - Adds the number of function evaluations spent
1320   in the line search to the running total.
1321 
1322    Input Parameters:
1323 +  tao - the TAO solver
1324 -  ls - the line search used in the optimization solver
1325 
1326    Level: developer
1327 
1328 .seealso: `TaoLineSearchApply()`
1329 @*/
1330 PetscErrorCode TaoAddLineSearchCounts(Tao tao)
1331 {
1332   PetscBool      flg;
1333   PetscInt       nfeval,ngeval,nfgeval;
1334 
1335   PetscFunctionBegin;
1336   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1337   if (tao->linesearch) {
1338     PetscCall(TaoLineSearchIsUsingTaoRoutines(tao->linesearch,&flg));
1339     if (!flg) {
1340       PetscCall(TaoLineSearchGetNumberFunctionEvaluations(tao->linesearch,&nfeval,&ngeval,&nfgeval));
1341       tao->nfuncs += nfeval;
1342       tao->ngrads += ngeval;
1343       tao->nfuncgrads += nfgeval;
1344     }
1345   }
1346   PetscFunctionReturn(0);
1347 }
1348 
1349 /*@
1350   TaoGetSolution - Returns the vector with the current TAO solution
1351 
1352   Not Collective
1353 
1354   Input Parameter:
1355 . tao - the Tao context
1356 
1357   Output Parameter:
1358 . X - the current solution
1359 
1360   Level: intermediate
1361 
1362   Note:  The returned vector will be the same object that was passed into TaoSetSolution()
1363 @*/
1364 PetscErrorCode TaoGetSolution(Tao tao, Vec *X)
1365 {
1366   PetscFunctionBegin;
1367   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1368   PetscValidPointer(X,2);
1369   *X = tao->solution;
1370   PetscFunctionReturn(0);
1371 }
1372 
1373 /*@
1374    TaoResetStatistics - Initialize the statistics used by TAO for all of the solvers.
1375    These statistics include the iteration number, residual norms, and convergence status.
1376    This routine gets called before solving each optimization problem.
1377 
1378    Collective on Tao
1379 
1380    Input Parameters:
1381 .  solver - the Tao context
1382 
1383    Level: developer
1384 
1385 .seealso: `TaoCreate()`, `TaoSolve()`
1386 @*/
1387 PetscErrorCode TaoResetStatistics(Tao tao)
1388 {
1389   PetscFunctionBegin;
1390   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1391   tao->niter        = 0;
1392   tao->nfuncs       = 0;
1393   tao->nfuncgrads   = 0;
1394   tao->ngrads       = 0;
1395   tao->nhess        = 0;
1396   tao->njac         = 0;
1397   tao->nconstraints = 0;
1398   tao->ksp_its      = 0;
1399   tao->ksp_tot_its  = 0;
1400   tao->reason       = TAO_CONTINUE_ITERATING;
1401   tao->residual     = 0.0;
1402   tao->cnorm        = 0.0;
1403   tao->step         = 0.0;
1404   tao->lsflag       = PETSC_FALSE;
1405   if (tao->hist_reset) tao->hist_len = 0;
1406   PetscFunctionReturn(0);
1407 }
1408 
1409 /*@C
1410   TaoSetUpdate - Sets the general-purpose update function called
1411   at the beginning of every iteration of the nonlinear solve. Specifically
1412   it is called at the top of every iteration, after the new solution and the gradient
1413   is determined, but before the Hessian is computed (if applicable).
1414 
1415   Logically Collective on Tao
1416 
1417   Input Parameters:
1418 + tao - The tao solver context
1419 - func - The function
1420 
1421   Calling sequence of func:
1422 $ func (Tao tao, PetscInt step);
1423 
1424 . step - The current step of the iteration
1425 
1426   Level: advanced
1427 
1428 .seealso `TaoSolve()`
1429 @*/
1430 PetscErrorCode TaoSetUpdate(Tao tao, PetscErrorCode (*func)(Tao, PetscInt, void*), void *ctx)
1431 {
1432   PetscFunctionBegin;
1433   PetscValidHeaderSpecific(tao, TAO_CLASSID,1);
1434   tao->ops->update = func;
1435   tao->user_update = ctx;
1436   PetscFunctionReturn(0);
1437 }
1438 
1439 /*@C
1440   TaoSetConvergenceTest - Sets the function that is to be used to test
1441   for convergence o fthe iterative minimization solution.  The new convergence
1442   testing routine will replace TAO's default convergence test.
1443 
1444   Logically Collective on Tao
1445 
1446   Input Parameters:
1447 + tao - the Tao object
1448 . conv - the routine to test for convergence
1449 - ctx - [optional] context for private data for the convergence routine
1450         (may be NULL)
1451 
1452   Calling sequence of conv:
1453 $   PetscErrorCode conv(Tao tao, void *ctx)
1454 
1455 + tao - the Tao object
1456 - ctx - [optional] convergence context
1457 
1458   Note: The new convergence testing routine should call TaoSetConvergedReason().
1459 
1460   Level: advanced
1461 
1462 .seealso: `TaoSetConvergedReason()`, `TaoGetSolutionStatus()`, `TaoGetTolerances()`, `TaoSetMonitor`
1463 
1464 @*/
1465 PetscErrorCode TaoSetConvergenceTest(Tao tao, PetscErrorCode (*conv)(Tao, void*), void *ctx)
1466 {
1467   PetscFunctionBegin;
1468   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1469   tao->ops->convergencetest = conv;
1470   tao->cnvP = ctx;
1471   PetscFunctionReturn(0);
1472 }
1473 
1474 /*@C
1475    TaoSetMonitor - Sets an ADDITIONAL function that is to be used at every
1476    iteration of the solver to display the iteration's
1477    progress.
1478 
1479    Logically Collective on Tao
1480 
1481    Input Parameters:
1482 +  tao - the Tao solver context
1483 .  mymonitor - monitoring routine
1484 -  mctx - [optional] user-defined context for private data for the
1485           monitor routine (may be NULL)
1486 
1487    Calling sequence of mymonitor:
1488 .vb
1489      PetscErrorCode mymonitor(Tao tao,void *mctx)
1490 .ve
1491 
1492 +    tao - the Tao solver context
1493 -    mctx - [optional] monitoring context
1494 
1495    Options Database Keys:
1496 +    -tao_monitor        - sets TaoMonitorDefault()
1497 .    -tao_smonitor       - sets short monitor
1498 .    -tao_cmonitor       - same as smonitor plus constraint norm
1499 .    -tao_view_solution   - view solution at each iteration
1500 .    -tao_view_gradient   - view gradient at each iteration
1501 .    -tao_view_ls_residual - view least-squares residual vector at each iteration
1502 -    -tao_cancelmonitors - cancels all monitors that have been hardwired into a code by calls to TaoSetMonitor(), but does not cancel those set via the options database.
1503 
1504    Notes:
1505    Several different monitoring routines may be set by calling
1506    TaoSetMonitor() multiple times; all will be called in the
1507    order in which they were set.
1508 
1509    Fortran Notes:
1510     Only one monitor function may be set
1511 
1512    Level: intermediate
1513 
1514 .seealso: `TaoMonitorDefault()`, `TaoCancelMonitors()`, `TaoSetDestroyRoutine()`
1515 @*/
1516 PetscErrorCode TaoSetMonitor(Tao tao, PetscErrorCode (*func)(Tao, void*), void *ctx,PetscErrorCode (*dest)(void**))
1517 {
1518   PetscInt       i;
1519   PetscBool      identical;
1520 
1521   PetscFunctionBegin;
1522   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1523   PetscCheck(tao->numbermonitors < MAXTAOMONITORS,PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Cannot attach another monitor -- max=%d",MAXTAOMONITORS);
1524 
1525   for (i=0; i<tao->numbermonitors;i++) {
1526     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))func,ctx,dest,(PetscErrorCode (*)(void))tao->monitor[i],tao->monitorcontext[i],tao->monitordestroy[i],&identical));
1527     if (identical) PetscFunctionReturn(0);
1528   }
1529   tao->monitor[tao->numbermonitors] = func;
1530   tao->monitorcontext[tao->numbermonitors] = (void*)ctx;
1531   tao->monitordestroy[tao->numbermonitors] = dest;
1532   ++tao->numbermonitors;
1533   PetscFunctionReturn(0);
1534 }
1535 
1536 /*@
1537    TaoCancelMonitors - Clears all the monitor functions for a Tao object.
1538 
1539    Logically Collective on Tao
1540 
1541    Input Parameters:
1542 .  tao - the Tao solver context
1543 
1544    Options Database:
1545 .  -tao_cancelmonitors - cancels all monitors that have been hardwired
1546     into a code by calls to TaoSetMonitor(), but does not cancel those
1547     set via the options database
1548 
1549    Notes:
1550    There is no way to clear one specific monitor from a Tao object.
1551 
1552    Level: advanced
1553 
1554 .seealso: `TaoMonitorDefault()`, `TaoSetMonitor()`
1555 @*/
1556 PetscErrorCode TaoCancelMonitors(Tao tao)
1557 {
1558   PetscInt       i;
1559 
1560   PetscFunctionBegin;
1561   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1562   for (i=0;i<tao->numbermonitors;i++) {
1563     if (tao->monitordestroy[i]) {
1564       PetscCall((*tao->monitordestroy[i])(&tao->monitorcontext[i]));
1565     }
1566   }
1567   tao->numbermonitors = 0;
1568   PetscFunctionReturn(0);
1569 }
1570 
1571 /*@
1572    TaoMonitorDefault - Default routine for monitoring progress of the
1573    Tao solvers (default).  This monitor prints the function value and gradient
1574    norm at each iteration.  It can be turned on from the command line using the
1575    -tao_monitor option
1576 
1577    Collective on Tao
1578 
1579    Input Parameters:
1580 +  tao - the Tao context
1581 -  ctx - PetscViewer context or NULL
1582 
1583    Options Database Keys:
1584 .  -tao_monitor - turn on default monitoring
1585 
1586    Level: advanced
1587 
1588 .seealso: `TaoDefaultSMonitor()`, `TaoSetMonitor()`
1589 @*/
1590 PetscErrorCode TaoMonitorDefault(Tao tao, void *ctx)
1591 {
1592   PetscInt       its, tabs;
1593   PetscReal      fct,gnorm;
1594   PetscViewer    viewer = (PetscViewer)ctx;
1595 
1596   PetscFunctionBegin;
1597   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1598   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1599   its = tao->niter;
1600   fct = tao->fc;
1601   gnorm = tao->residual;
1602   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1603   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1604   if (its == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
1605      PetscCall(PetscViewerASCIIPrintf(viewer,"  Iteration information for %s solve.\n",((PetscObject)tao)->prefix));
1606      tao->header_printed = PETSC_TRUE;
1607    }
1608   PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " TAO,",its));
1609   PetscCall(PetscViewerASCIIPrintf(viewer,"  Function value: %g,",(double)fct));
1610   if (gnorm >= PETSC_INFINITY) {
1611     PetscCall(PetscViewerASCIIPrintf(viewer,"  Residual: Inf \n"));
1612   } else {
1613     PetscCall(PetscViewerASCIIPrintf(viewer,"  Residual: %g \n",(double)gnorm));
1614   }
1615   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1616   PetscFunctionReturn(0);
1617 }
1618 
1619 /*@
1620    TaoDefaultGMonitor - Default routine for monitoring progress of the
1621    Tao solvers (default) with extra detail on the globalization method.
1622    This monitor prints the function value and gradient norm at each
1623    iteration, as well as the step size and trust radius. Note that the
1624    step size and trust radius may be the same for some algorithms.
1625    It can be turned on from the command line using the
1626    -tao_gmonitor option
1627 
1628    Collective on Tao
1629 
1630    Input Parameters:
1631 +  tao - the Tao context
1632 -  ctx - PetscViewer context or NULL
1633 
1634    Options Database Keys:
1635 .  -tao_gmonitor - turn on monitoring with globalization information
1636 
1637    Level: advanced
1638 
1639 .seealso: `TaoDefaultSMonitor()`, `TaoSetMonitor()`
1640 @*/
1641 PetscErrorCode TaoDefaultGMonitor(Tao tao, void *ctx)
1642 {
1643   PetscInt       its, tabs;
1644   PetscReal      fct,gnorm,stp,tr;
1645   PetscViewer    viewer = (PetscViewer)ctx;
1646 
1647   PetscFunctionBegin;
1648   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1649   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1650   its = tao->niter;
1651   fct = tao->fc;
1652   gnorm = tao->residual;
1653   stp = tao->step;
1654   tr = tao->trust;
1655   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1656   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1657   if (its == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
1658      PetscCall(PetscViewerASCIIPrintf(viewer,"  Iteration information for %s solve.\n",((PetscObject)tao)->prefix));
1659      tao->header_printed = PETSC_TRUE;
1660    }
1661   PetscCall(PetscViewerASCIIPrintf(viewer,"%3" PetscInt_FMT " TAO,",its));
1662   PetscCall(PetscViewerASCIIPrintf(viewer,"  Function value: %g,",(double)fct));
1663   if (gnorm >= PETSC_INFINITY) {
1664     PetscCall(PetscViewerASCIIPrintf(viewer,"  Residual: Inf,"));
1665   } else {
1666     PetscCall(PetscViewerASCIIPrintf(viewer,"  Residual: %g,",(double)gnorm));
1667   }
1668   PetscCall(PetscViewerASCIIPrintf(viewer,"  Step: %g,  Trust: %g\n",(double)stp,(double)tr));
1669   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1670   PetscFunctionReturn(0);
1671 }
1672 
1673 /*@
1674    TaoDefaultSMonitor - Default routine for monitoring progress of the
1675    solver. Same as TaoMonitorDefault() except
1676    it prints fewer digits of the residual as the residual gets smaller.
1677    This is because the later digits are meaningless and are often
1678    different on different machines; by using this routine different
1679    machines will usually generate the same output. It can be turned on
1680    by using the -tao_smonitor option
1681 
1682    Collective on Tao
1683 
1684    Input Parameters:
1685 +  tao - the Tao context
1686 -  ctx - PetscViewer context of type ASCII
1687 
1688    Options Database Keys:
1689 .  -tao_smonitor - turn on default short monitoring
1690 
1691    Level: advanced
1692 
1693 .seealso: `TaoMonitorDefault()`, `TaoSetMonitor()`
1694 @*/
1695 PetscErrorCode TaoDefaultSMonitor(Tao tao, void *ctx)
1696 {
1697   PetscInt       its, tabs;
1698   PetscReal      fct,gnorm;
1699   PetscViewer    viewer = (PetscViewer)ctx;
1700 
1701   PetscFunctionBegin;
1702   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1703   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1704   its = tao->niter;
1705   fct = tao->fc;
1706   gnorm = tao->residual;
1707   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1708   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1709   PetscCall(PetscViewerASCIIPrintf(viewer,"iter = %3" PetscInt_FMT ",",its));
1710   PetscCall(PetscViewerASCIIPrintf(viewer," Function value %g,",(double)fct));
1711   if (gnorm >= PETSC_INFINITY) {
1712     PetscCall(PetscViewerASCIIPrintf(viewer," Residual: Inf \n"));
1713   } else if (gnorm > 1.e-6) {
1714     PetscCall(PetscViewerASCIIPrintf(viewer," Residual: %g \n",(double)gnorm));
1715   } else if (gnorm > 1.e-11) {
1716     PetscCall(PetscViewerASCIIPrintf(viewer," Residual: < 1.0e-6 \n"));
1717   } else {
1718     PetscCall(PetscViewerASCIIPrintf(viewer," Residual: < 1.0e-11 \n"));
1719   }
1720   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1721   PetscFunctionReturn(0);
1722 }
1723 
1724 /*@
1725    TaoDefaultCMonitor - same as TaoMonitorDefault() except
1726    it prints the norm of the constraints function. It can be turned on
1727    from the command line using the -tao_cmonitor option
1728 
1729    Collective on Tao
1730 
1731    Input Parameters:
1732 +  tao - the Tao context
1733 -  ctx - PetscViewer context or NULL
1734 
1735    Options Database Keys:
1736 .  -tao_cmonitor - monitor the constraints
1737 
1738    Level: advanced
1739 
1740 .seealso: `TaoMonitorDefault()`, `TaoSetMonitor()`
1741 @*/
1742 PetscErrorCode TaoDefaultCMonitor(Tao tao, void *ctx)
1743 {
1744   PetscInt       its, tabs;
1745   PetscReal      fct,gnorm;
1746   PetscViewer    viewer = (PetscViewer)ctx;
1747 
1748   PetscFunctionBegin;
1749   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1750   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1751   its = tao->niter;
1752   fct = tao->fc;
1753   gnorm = tao->residual;
1754   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1755   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1756   PetscCall(PetscViewerASCIIPrintf(viewer,"iter = %" PetscInt_FMT ",",its));
1757   PetscCall(PetscViewerASCIIPrintf(viewer," Function value: %g,",(double)fct));
1758   PetscCall(PetscViewerASCIIPrintf(viewer,"  Residual: %g ",(double)gnorm));
1759   PetscCall(PetscViewerASCIIPrintf(viewer,"  Constraint: %g \n",(double)tao->cnorm));
1760   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1761   PetscFunctionReturn(0);
1762 }
1763 
1764 /*@C
1765    TaoSolutionMonitor - Views the solution at each iteration
1766    It can be turned on from the command line using the
1767    -tao_view_solution option
1768 
1769    Collective on Tao
1770 
1771    Input Parameters:
1772 +  tao - the Tao context
1773 -  ctx - PetscViewer context or NULL
1774 
1775    Options Database Keys:
1776 .  -tao_view_solution - view the solution
1777 
1778    Level: advanced
1779 
1780 .seealso: `TaoDefaultSMonitor()`, `TaoSetMonitor()`
1781 @*/
1782 PetscErrorCode TaoSolutionMonitor(Tao tao, void *ctx)
1783 {
1784   PetscViewer    viewer  = (PetscViewer)ctx;
1785 
1786   PetscFunctionBegin;
1787   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1788   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1789   PetscCall(VecView(tao->solution,viewer));
1790   PetscFunctionReturn(0);
1791 }
1792 
1793 /*@C
1794    TaoGradientMonitor - Views the gradient at each iteration
1795    It can be turned on from the command line using the
1796    -tao_view_gradient option
1797 
1798    Collective on Tao
1799 
1800    Input Parameters:
1801 +  tao - the Tao context
1802 -  ctx - PetscViewer context or NULL
1803 
1804    Options Database Keys:
1805 .  -tao_view_gradient - view the gradient at each iteration
1806 
1807    Level: advanced
1808 
1809 .seealso: `TaoDefaultSMonitor()`, `TaoSetMonitor()`
1810 @*/
1811 PetscErrorCode TaoGradientMonitor(Tao tao, void *ctx)
1812 {
1813   PetscViewer    viewer = (PetscViewer)ctx;
1814 
1815   PetscFunctionBegin;
1816   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1817   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1818   PetscCall(VecView(tao->gradient,viewer));
1819   PetscFunctionReturn(0);
1820 }
1821 
1822 /*@C
1823    TaoStepDirectionMonitor - Views the step-direction at each iteration
1824 
1825    Collective on Tao
1826 
1827    Input Parameters:
1828 +  tao - the Tao context
1829 -  ctx - PetscViewer context or NULL
1830 
1831    Options Database Keys:
1832 .  -tao_view_gradient - view the gradient at each iteration
1833 
1834    Level: advanced
1835 
1836 .seealso: `TaoDefaultSMonitor()`, `TaoSetMonitor()`
1837 @*/
1838 PetscErrorCode TaoStepDirectionMonitor(Tao tao, void *ctx)
1839 {
1840   PetscViewer    viewer = (PetscViewer)ctx;
1841 
1842   PetscFunctionBegin;
1843   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1844   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1845   PetscCall(VecView(tao->stepdirection,viewer));
1846   PetscFunctionReturn(0);
1847 }
1848 
1849 /*@C
1850    TaoDrawSolutionMonitor - Plots the solution at each iteration
1851    It can be turned on from the command line using the
1852    -tao_draw_solution option
1853 
1854    Collective on Tao
1855 
1856    Input Parameters:
1857 +  tao - the Tao context
1858 -  ctx - TaoMonitorDraw context
1859 
1860    Options Database Keys:
1861 .  -tao_draw_solution - draw the solution at each iteration
1862 
1863    Level: advanced
1864 
1865 .seealso: `TaoSolutionMonitor()`, `TaoSetMonitor()`, `TaoDrawGradientMonitor`
1866 @*/
1867 PetscErrorCode TaoDrawSolutionMonitor(Tao tao, void *ctx)
1868 {
1869   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;
1870 
1871   PetscFunctionBegin;
1872   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1873   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(0);
1874   PetscCall(VecView(tao->solution,ictx->viewer));
1875   PetscFunctionReturn(0);
1876 }
1877 
1878 /*@C
1879    TaoDrawGradientMonitor - Plots the gradient at each iteration
1880    It can be turned on from the command line using the
1881    -tao_draw_gradient option
1882 
1883    Collective on Tao
1884 
1885    Input Parameters:
1886 +  tao - the Tao context
1887 -  ctx - PetscViewer context
1888 
1889    Options Database Keys:
1890 .  -tao_draw_gradient - draw the gradient at each iteration
1891 
1892    Level: advanced
1893 
1894 .seealso: `TaoGradientMonitor()`, `TaoSetMonitor()`, `TaoDrawSolutionMonitor`
1895 @*/
1896 PetscErrorCode TaoDrawGradientMonitor(Tao tao, void *ctx)
1897 {
1898   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;
1899 
1900   PetscFunctionBegin;
1901   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1902   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(0);
1903   PetscCall(VecView(tao->gradient,ictx->viewer));
1904   PetscFunctionReturn(0);
1905 }
1906 
1907 /*@C
1908    TaoDrawStepMonitor - Plots the step direction at each iteration
1909 
1910    Collective on Tao
1911 
1912    Input Parameters:
1913 +  tao - the Tao context
1914 -  ctx - PetscViewer context
1915 
1916    Options Database Keys:
1917 .  -tao_draw_step - draw the step direction at each iteration
1918 
1919    Level: advanced
1920 
1921 .seealso: `TaoSetMonitor()`, `TaoDrawSolutionMonitor`
1922 @*/
1923 PetscErrorCode TaoDrawStepMonitor(Tao tao, void *ctx)
1924 {
1925   PetscViewer    viewer = (PetscViewer)ctx;
1926 
1927   PetscFunctionBegin;
1928   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1929   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1930   PetscCall(VecView(tao->stepdirection,viewer));
1931   PetscFunctionReturn(0);
1932 }
1933 
1934 /*@C
1935    TaoResidualMonitor - Views the least-squares residual at each iteration
1936 
1937    Collective on Tao
1938 
1939    Input Parameters:
1940 +  tao - the Tao context
1941 -  ctx - PetscViewer context or NULL
1942 
1943    Options Database Keys:
1944 .  -tao_view_ls_residual - view the least-squares residual at each iteration
1945 
1946    Level: advanced
1947 
1948 .seealso: `TaoDefaultSMonitor()`, `TaoSetMonitor()`
1949 @*/
1950 PetscErrorCode TaoResidualMonitor(Tao tao, void *ctx)
1951 {
1952   PetscViewer    viewer  = (PetscViewer)ctx;
1953 
1954   PetscFunctionBegin;
1955   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1956   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1957   PetscCall(VecView(tao->ls_res,viewer));
1958   PetscFunctionReturn(0);
1959 }
1960 
1961 /*@
1962    TaoDefaultConvergenceTest - Determines whether the solver should continue iterating
1963    or terminate.
1964 
1965    Collective on Tao
1966 
1967    Input Parameters:
1968 +  tao - the Tao context
1969 -  dummy - unused dummy context
1970 
1971    Output Parameter:
1972 .  reason - for terminating
1973 
1974    Notes:
1975    This routine checks the residual in the optimality conditions, the
1976    relative residual in the optimity conditions, the number of function
1977    evaluations, and the function value to test convergence.  Some
1978    solvers may use different convergence routines.
1979 
1980    Level: developer
1981 
1982 .seealso: `TaoSetTolerances()`, `TaoGetConvergedReason()`, `TaoSetConvergedReason()`
1983 @*/
1984 
1985 PetscErrorCode TaoDefaultConvergenceTest(Tao tao,void *dummy)
1986 {
1987   PetscInt           niter=tao->niter, nfuncs=PetscMax(tao->nfuncs,tao->nfuncgrads);
1988   PetscInt           max_funcs=tao->max_funcs;
1989   PetscReal          gnorm=tao->residual, gnorm0=tao->gnorm0;
1990   PetscReal          f=tao->fc, steptol=tao->steptol,trradius=tao->step;
1991   PetscReal          gatol=tao->gatol,grtol=tao->grtol,gttol=tao->gttol;
1992   PetscReal          catol=tao->catol,crtol=tao->crtol;
1993   PetscReal          fmin=tao->fmin, cnorm=tao->cnorm;
1994   TaoConvergedReason reason=tao->reason;
1995 
1996   PetscFunctionBegin;
1997   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1998   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
1999 
2000   if (PetscIsInfOrNanReal(f)) {
2001     PetscCall(PetscInfo(tao,"Failed to converged, function value is Inf or NaN\n"));
2002     reason = TAO_DIVERGED_NAN;
2003   } else if (f <= fmin && cnorm <=catol) {
2004     PetscCall(PetscInfo(tao,"Converged due to function value %g < minimum function value %g\n", (double)f,(double)fmin));
2005     reason = TAO_CONVERGED_MINF;
2006   } else if (gnorm<= gatol && cnorm <=catol) {
2007     PetscCall(PetscInfo(tao,"Converged due to residual norm ||g(X)||=%g < %g\n",(double)gnorm,(double)gatol));
2008     reason = TAO_CONVERGED_GATOL;
2009   } else if (f!=0 && PetscAbsReal(gnorm/f) <= grtol && cnorm <= crtol) {
2010     PetscCall(PetscInfo(tao,"Converged due to residual ||g(X)||/|f(X)| =%g < %g\n",(double)(gnorm/f),(double)grtol));
2011     reason = TAO_CONVERGED_GRTOL;
2012   } else if (gnorm0 != 0 && ((gttol == 0 && gnorm == 0) || gnorm/gnorm0 < gttol) && cnorm <= crtol) {
2013     PetscCall(PetscInfo(tao,"Converged due to relative residual norm ||g(X)||/||g(X0)|| = %g < %g\n",(double)(gnorm/gnorm0),(double)gttol));
2014     reason = TAO_CONVERGED_GTTOL;
2015   } else if (max_funcs >=0 && nfuncs > max_funcs) {
2016     PetscCall(PetscInfo(tao,"Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", nfuncs,max_funcs));
2017     reason = TAO_DIVERGED_MAXFCN;
2018   } else if (tao->lsflag != 0) {
2019     PetscCall(PetscInfo(tao,"Tao Line Search failure.\n"));
2020     reason = TAO_DIVERGED_LS_FAILURE;
2021   } else if (trradius < steptol && niter > 0) {
2022     PetscCall(PetscInfo(tao,"Trust region/step size too small: %g < %g\n", (double)trradius,(double)steptol));
2023     reason = TAO_CONVERGED_STEPTOL;
2024   } else if (niter >= tao->max_it) {
2025     PetscCall(PetscInfo(tao,"Exceeded maximum number of iterations: %" PetscInt_FMT " > %" PetscInt_FMT "\n",niter,tao->max_it));
2026     reason = TAO_DIVERGED_MAXITS;
2027   } else {
2028     reason = TAO_CONTINUE_ITERATING;
2029   }
2030   tao->reason = reason;
2031   PetscFunctionReturn(0);
2032 }
2033 
2034 /*@C
2035    TaoSetOptionsPrefix - Sets the prefix used for searching for all
2036    TAO options in the database.
2037 
2038    Logically Collective on Tao
2039 
2040    Input Parameters:
2041 +  tao - the Tao context
2042 -  prefix - the prefix string to prepend to all TAO option requests
2043 
2044    Notes:
2045    A hyphen (-) must NOT be given at the beginning of the prefix name.
2046    The first character of all runtime options is AUTOMATICALLY the hyphen.
2047 
2048    For example, to distinguish between the runtime options for two
2049    different TAO solvers, one could call
2050 .vb
2051       TaoSetOptionsPrefix(tao1,"sys1_")
2052       TaoSetOptionsPrefix(tao2,"sys2_")
2053 .ve
2054 
2055    This would enable use of different options for each system, such as
2056 .vb
2057       -sys1_tao_method blmvm -sys1_tao_grtol 1.e-3
2058       -sys2_tao_method lmvm  -sys2_tao_grtol 1.e-4
2059 .ve
2060 
2061    Level: advanced
2062 
2063 .seealso: `TaoAppendOptionsPrefix()`, `TaoGetOptionsPrefix()`
2064 @*/
2065 
2066 PetscErrorCode TaoSetOptionsPrefix(Tao tao, const char p[])
2067 {
2068   PetscFunctionBegin;
2069   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2070   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tao,p));
2071   if (tao->linesearch) {
2072     PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,p));
2073   }
2074   if (tao->ksp) {
2075     PetscCall(KSPSetOptionsPrefix(tao->ksp,p));
2076   }
2077   PetscFunctionReturn(0);
2078 }
2079 
2080 /*@C
2081    TaoAppendOptionsPrefix - Appends to the prefix used for searching for all
2082    TAO options in the database.
2083 
2084    Logically Collective on Tao
2085 
2086    Input Parameters:
2087 +  tao - the Tao solver context
2088 -  prefix - the prefix string to prepend to all TAO option requests
2089 
2090    Notes:
2091    A hyphen (-) must NOT be given at the beginning of the prefix name.
2092    The first character of all runtime options is AUTOMATICALLY the hyphen.
2093 
2094    Level: advanced
2095 
2096 .seealso: `TaoSetOptionsPrefix()`, `TaoGetOptionsPrefix()`
2097 @*/
2098 PetscErrorCode TaoAppendOptionsPrefix(Tao tao, const char p[])
2099 {
2100   PetscFunctionBegin;
2101   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2102   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao,p));
2103   if (tao->linesearch) {
2104     PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao->linesearch,p));
2105   }
2106   if (tao->ksp) {
2107     PetscCall(KSPAppendOptionsPrefix(tao->ksp,p));
2108   }
2109   PetscFunctionReturn(0);
2110 }
2111 
2112 /*@C
2113   TaoGetOptionsPrefix - Gets the prefix used for searching for all
2114   TAO options in the database
2115 
2116   Not Collective
2117 
2118   Input Parameters:
2119 . tao - the Tao context
2120 
2121   Output Parameters:
2122 . prefix - pointer to the prefix string used is returned
2123 
2124   Notes:
2125     On the fortran side, the user should pass in a string 'prefix' of
2126   sufficient length to hold the prefix.
2127 
2128   Level: advanced
2129 
2130 .seealso: `TaoSetOptionsPrefix()`, `TaoAppendOptionsPrefix()`
2131 @*/
2132 PetscErrorCode TaoGetOptionsPrefix(Tao tao, const char *p[])
2133 {
2134   PetscFunctionBegin;
2135   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2136   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)tao,p));
2137   PetscFunctionReturn(0);
2138 }
2139 
2140 /*@C
2141    TaoSetType - Sets the method for the unconstrained minimization solver.
2142 
2143    Collective on Tao
2144 
2145    Input Parameters:
2146 +  solver - the Tao solver context
2147 -  type - a known method
2148 
2149    Options Database Key:
2150 .  -tao_type <type> - Sets the method; use -help for a list
2151    of available methods (for instance, "-tao_type lmvm" or "-tao_type tron")
2152 
2153    Available methods include:
2154 +    nls - Newton's method with line search for unconstrained minimization
2155 .    ntr - Newton's method with trust region for unconstrained minimization
2156 .    ntl - Newton's method with trust region, line search for unconstrained minimization
2157 .    lmvm - Limited memory variable metric method for unconstrained minimization
2158 .    cg - Nonlinear conjugate gradient method for unconstrained minimization
2159 .    nm - Nelder-Mead algorithm for derivate-free unconstrained minimization
2160 .    tron - Newton Trust Region method for bound constrained minimization
2161 .    gpcg - Newton Trust Region method for quadratic bound constrained minimization
2162 .    blmvm - Limited memory variable metric method for bound constrained minimization
2163 -    pounders - Model-based algorithm pounder extended for nonlinear least squares
2164 
2165   Level: intermediate
2166 
2167 .seealso: `TaoCreate()`, `TaoGetType()`, `TaoType`
2168 
2169 @*/
2170 PetscErrorCode TaoSetType(Tao tao, TaoType type)
2171 {
2172   PetscErrorCode (*create_xxx)(Tao);
2173   PetscBool      issame;
2174 
2175   PetscFunctionBegin;
2176   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2177 
2178   PetscCall(PetscObjectTypeCompare((PetscObject)tao,type,&issame));
2179   if (issame) PetscFunctionReturn(0);
2180 
2181   PetscCall(PetscFunctionListFind(TaoList, type, (void(**)(void))&create_xxx));
2182   PetscCheck(create_xxx,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Tao type %s",type);
2183 
2184   /* Destroy the existing solver information */
2185   if (tao->ops->destroy) PetscCall((*tao->ops->destroy)(tao));
2186   PetscCall(TaoLineSearchDestroy(&tao->linesearch));
2187   tao->ops->setup          = NULL;
2188   tao->ops->solve          = NULL;
2189   tao->ops->view           = NULL;
2190   tao->ops->setfromoptions = NULL;
2191   tao->ops->destroy        = NULL;
2192 
2193   tao->setupcalled = PETSC_FALSE;
2194 
2195   PetscCall((*create_xxx)(tao));
2196   PetscCall(PetscObjectChangeTypeName((PetscObject)tao,type));
2197   PetscFunctionReturn(0);
2198 }
2199 
2200 /*MC
2201    TaoRegister - Adds a method to the TAO package for unconstrained minimization.
2202 
2203    Synopsis:
2204    TaoRegister(char *name_solver,char *path,char *name_Create,PetscErrorCode (*routine_Create)(Tao))
2205 
2206    Not collective
2207 
2208    Input Parameters:
2209 +  sname - name of a new user-defined solver
2210 -  func - routine to Create method context
2211 
2212    Notes:
2213    TaoRegister() may be called multiple times to add several user-defined solvers.
2214 
2215    Sample usage:
2216 .vb
2217    TaoRegister("my_solver",MySolverCreate);
2218 .ve
2219 
2220    Then, your solver can be chosen with the procedural interface via
2221 $     TaoSetType(tao,"my_solver")
2222    or at runtime via the option
2223 $     -tao_type my_solver
2224 
2225    Level: advanced
2226 
2227 .seealso: `TaoRegisterAll()`, `TaoRegisterDestroy()`
2228 M*/
2229 PetscErrorCode TaoRegister(const char sname[], PetscErrorCode (*func)(Tao))
2230 {
2231   PetscFunctionBegin;
2232   PetscCall(TaoInitializePackage());
2233   PetscCall(PetscFunctionListAdd(&TaoList,sname,(void (*)(void))func));
2234   PetscFunctionReturn(0);
2235 }
2236 
2237 /*@C
2238    TaoRegisterDestroy - Frees the list of minimization solvers that were
2239    registered by TaoRegisterDynamic().
2240 
2241    Not Collective
2242 
2243    Level: advanced
2244 
2245 .seealso: `TaoRegisterAll()`, `TaoRegister()`
2246 @*/
2247 PetscErrorCode TaoRegisterDestroy(void)
2248 {
2249   PetscFunctionBegin;
2250   PetscCall(PetscFunctionListDestroy(&TaoList));
2251   TaoRegisterAllCalled = PETSC_FALSE;
2252   PetscFunctionReturn(0);
2253 }
2254 
2255 /*@
2256    TaoGetIterationNumber - Gets the number of Tao iterations completed
2257    at this time.
2258 
2259    Not Collective
2260 
2261    Input Parameter:
2262 .  tao - Tao context
2263 
2264    Output Parameter:
2265 .  iter - iteration number
2266 
2267    Notes:
2268    For example, during the computation of iteration 2 this would return 1.
2269 
2270    Level: intermediate
2271 
2272 .seealso: `TaoGetLinearSolveIterations()`, `TaoGetResidualNorm()`, `TaoGetObjective()`
2273 @*/
2274 PetscErrorCode  TaoGetIterationNumber(Tao tao,PetscInt *iter)
2275 {
2276   PetscFunctionBegin;
2277   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2278   PetscValidIntPointer(iter,2);
2279   *iter = tao->niter;
2280   PetscFunctionReturn(0);
2281 }
2282 
2283 /*@
2284    TaoGetResidualNorm - Gets the current value of the norm of the residual
2285    at this time.
2286 
2287    Not Collective
2288 
2289    Input Parameter:
2290 .  tao - Tao context
2291 
2292    Output Parameter:
2293 .  value - the current value
2294 
2295    Level: intermediate
2296 
2297    Developer Note: This is the 2-norm of the residual, we cannot use TaoGetGradientNorm() because that has
2298                    a different meaning. For some reason Tao sometimes calls the gradient the residual.
2299 
2300 .seealso: `TaoGetLinearSolveIterations()`, `TaoGetIterationNumber()`, `TaoGetObjective()`
2301 @*/
2302 PetscErrorCode TaoGetResidualNorm(Tao tao,PetscReal *value)
2303 {
2304   PetscFunctionBegin;
2305   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2306   PetscValidRealPointer(value,2);
2307   *value = tao->residual;
2308   PetscFunctionReturn(0);
2309 }
2310 
2311 /*@
2312    TaoSetIterationNumber - Sets the current iteration number.
2313 
2314    Logically Collective on Tao
2315 
2316    Input Parameters:
2317 +  tao - Tao context
2318 -  iter - iteration number
2319 
2320    Level: developer
2321 
2322 .seealso: `TaoGetLinearSolveIterations()`
2323 @*/
2324 PetscErrorCode  TaoSetIterationNumber(Tao tao,PetscInt iter)
2325 {
2326   PetscFunctionBegin;
2327   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2328   PetscValidLogicalCollectiveInt(tao,iter,2);
2329   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2330   tao->niter = iter;
2331   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2332   PetscFunctionReturn(0);
2333 }
2334 
2335 /*@
2336    TaoGetTotalIterationNumber - Gets the total number of Tao iterations
2337    completed. This number keeps accumulating if multiple solves
2338    are called with the Tao object.
2339 
2340    Not Collective
2341 
2342    Input Parameter:
2343 .  tao - Tao context
2344 
2345    Output Parameter:
2346 .  iter - iteration number
2347 
2348    Notes:
2349    The total iteration count is updated after each solve, if there is a current
2350    TaoSolve() in progress then those iterations are not yet counted.
2351 
2352    Level: intermediate
2353 
2354 .seealso: `TaoGetLinearSolveIterations()`
2355 @*/
2356 PetscErrorCode  TaoGetTotalIterationNumber(Tao tao,PetscInt *iter)
2357 {
2358   PetscFunctionBegin;
2359   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2360   PetscValidIntPointer(iter,2);
2361   *iter = tao->ntotalits;
2362   PetscFunctionReturn(0);
2363 }
2364 
2365 /*@
2366    TaoSetTotalIterationNumber - Sets the current total iteration number.
2367 
2368    Logically Collective on Tao
2369 
2370    Input Parameters:
2371 +  tao - Tao context
2372 -  iter - iteration number
2373 
2374    Level: developer
2375 
2376 .seealso: `TaoGetLinearSolveIterations()`
2377 @*/
2378 PetscErrorCode  TaoSetTotalIterationNumber(Tao tao,PetscInt iter)
2379 {
2380   PetscFunctionBegin;
2381   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2382   PetscValidLogicalCollectiveInt(tao,iter,2);
2383   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2384   tao->ntotalits = iter;
2385   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2386   PetscFunctionReturn(0);
2387 }
2388 
2389 /*@
2390   TaoSetConvergedReason - Sets the termination flag on a Tao object
2391 
2392   Logically Collective on Tao
2393 
2394   Input Parameters:
2395 + tao - the Tao context
2396 - reason - one of
2397 $     TAO_CONVERGED_ATOL (2),
2398 $     TAO_CONVERGED_RTOL (3),
2399 $     TAO_CONVERGED_STEPTOL (4),
2400 $     TAO_CONVERGED_MINF (5),
2401 $     TAO_CONVERGED_USER (6),
2402 $     TAO_DIVERGED_MAXITS (-2),
2403 $     TAO_DIVERGED_NAN (-4),
2404 $     TAO_DIVERGED_MAXFCN (-5),
2405 $     TAO_DIVERGED_LS_FAILURE (-6),
2406 $     TAO_DIVERGED_TR_REDUCTION (-7),
2407 $     TAO_DIVERGED_USER (-8),
2408 $     TAO_CONTINUE_ITERATING (0)
2409 
2410    Level: intermediate
2411 
2412 @*/
2413 PetscErrorCode TaoSetConvergedReason(Tao tao, TaoConvergedReason reason)
2414 {
2415   PetscFunctionBegin;
2416   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2417   PetscValidLogicalCollectiveEnum(tao,reason,2);
2418   tao->reason = reason;
2419   PetscFunctionReturn(0);
2420 }
2421 
2422 /*@
2423    TaoGetConvergedReason - Gets the reason the Tao iteration was stopped.
2424 
2425    Not Collective
2426 
2427    Input Parameter:
2428 .  tao - the Tao solver context
2429 
2430    Output Parameter:
2431 .  reason - one of
2432 $  TAO_CONVERGED_GATOL (3)           ||g(X)|| < gatol
2433 $  TAO_CONVERGED_GRTOL (4)           ||g(X)|| / f(X)  < grtol
2434 $  TAO_CONVERGED_GTTOL (5)           ||g(X)|| / ||g(X0)|| < gttol
2435 $  TAO_CONVERGED_STEPTOL (6)         step size small
2436 $  TAO_CONVERGED_MINF (7)            F < F_min
2437 $  TAO_CONVERGED_USER (8)            User defined
2438 $  TAO_DIVERGED_MAXITS (-2)          its > maxits
2439 $  TAO_DIVERGED_NAN (-4)             Numerical problems
2440 $  TAO_DIVERGED_MAXFCN (-5)          fevals > max_funcsals
2441 $  TAO_DIVERGED_LS_FAILURE (-6)      line search failure
2442 $  TAO_DIVERGED_TR_REDUCTION (-7)    trust region failure
2443 $  TAO_DIVERGED_USER(-8)             (user defined)
2444 $  TAO_CONTINUE_ITERATING (0)
2445 
2446    where
2447 +  X - current solution
2448 .  X0 - initial guess
2449 .  f(X) - current function value
2450 .  f(X*) - true solution (estimated)
2451 .  g(X) - current gradient
2452 .  its - current iterate number
2453 .  maxits - maximum number of iterates
2454 .  fevals - number of function evaluations
2455 -  max_funcsals - maximum number of function evaluations
2456 
2457    Level: intermediate
2458 
2459 .seealso: `TaoSetConvergenceTest()`, `TaoSetTolerances()`
2460 
2461 @*/
2462 PetscErrorCode TaoGetConvergedReason(Tao tao, TaoConvergedReason *reason)
2463 {
2464   PetscFunctionBegin;
2465   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2466   PetscValidPointer(reason,2);
2467   *reason = tao->reason;
2468   PetscFunctionReturn(0);
2469 }
2470 
2471 /*@
2472    TaoGetSolutionStatus - Get the current iterate, objective value,
2473    residual, infeasibility, and termination
2474 
2475    Not Collective
2476 
2477    Input Parameter:
2478 .  tao - the Tao context
2479 
2480    Output Parameters:
2481 +  iterate - the current iterate number (>=0)
2482 .  f - the current function value
2483 .  gnorm - the square of the gradient norm, duality gap, or other measure indicating distance from optimality.
2484 .  cnorm - the infeasibility of the current solution with regard to the constraints.
2485 .  xdiff - the step length or trust region radius of the most recent iterate.
2486 -  reason - The termination reason, which can equal TAO_CONTINUE_ITERATING
2487 
2488    Level: intermediate
2489 
2490    Note:
2491    TAO returns the values set by the solvers in the routine TaoMonitor().
2492 
2493    Note:
2494    If any of the output arguments are set to NULL, no corresponding value will be returned.
2495 
2496 .seealso: `TaoMonitor()`, `TaoGetConvergedReason()`
2497 @*/
2498 PetscErrorCode TaoGetSolutionStatus(Tao tao, PetscInt *its, PetscReal *f, PetscReal *gnorm, PetscReal *cnorm, PetscReal *xdiff, TaoConvergedReason *reason)
2499 {
2500   PetscFunctionBegin;
2501   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2502   if (its) *its = tao->niter;
2503   if (f) *f = tao->fc;
2504   if (gnorm) *gnorm = tao->residual;
2505   if (cnorm) *cnorm = tao->cnorm;
2506   if (reason) *reason = tao->reason;
2507   if (xdiff) *xdiff = tao->step;
2508   PetscFunctionReturn(0);
2509 }
2510 
2511 /*@C
2512    TaoGetType - Gets the current Tao algorithm.
2513 
2514    Not Collective
2515 
2516    Input Parameter:
2517 .  tao - the Tao solver context
2518 
2519    Output Parameter:
2520 .  type - Tao method
2521 
2522    Level: intermediate
2523 
2524 @*/
2525 PetscErrorCode TaoGetType(Tao tao,TaoType *type)
2526 {
2527   PetscFunctionBegin;
2528   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2529   PetscValidPointer(type,2);
2530   *type = ((PetscObject)tao)->type_name;
2531   PetscFunctionReturn(0);
2532 }
2533 
2534 /*@C
2535   TaoMonitor - Monitor the solver and the current solution.  This
2536   routine will record the iteration number and residual statistics,
2537   call any monitors specified by the user, and calls the convergence-check routine.
2538 
2539    Input Parameters:
2540 +  tao - the Tao context
2541 .  its - the current iterate number (>=0)
2542 .  f - the current objective function value
2543 .  res - the gradient norm, square root of the duality gap, or other measure indicating distince from optimality.  This measure will be recorded and
2544           used for some termination tests.
2545 .  cnorm - the infeasibility of the current solution with regard to the constraints.
2546 -  steplength - multiple of the step direction added to the previous iterate.
2547 
2548    Output Parameters:
2549 .  reason - The termination reason, which can equal TAO_CONTINUE_ITERATING
2550 
2551    Options Database Key:
2552 .  -tao_monitor - Use the default monitor, which prints statistics to standard output
2553 
2554 .seealso `TaoGetConvergedReason()`, `TaoMonitorDefault()`, `TaoSetMonitor()`
2555 
2556    Level: developer
2557 
2558 @*/
2559 PetscErrorCode TaoMonitor(Tao tao, PetscInt its, PetscReal f, PetscReal res, PetscReal cnorm, PetscReal steplength)
2560 {
2561   PetscInt       i;
2562 
2563   PetscFunctionBegin;
2564   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2565   tao->fc = f;
2566   tao->residual = res;
2567   tao->cnorm = cnorm;
2568   tao->step = steplength;
2569   if (!its) {
2570     tao->cnorm0 = cnorm;
2571     tao->gnorm0 = res;
2572   }
2573   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(res),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
2574   for (i=0;i<tao->numbermonitors;i++) {
2575     PetscCall((*tao->monitor[i])(tao,tao->monitorcontext[i]));
2576   }
2577   PetscFunctionReturn(0);
2578 }
2579 
2580 /*@
2581    TaoSetConvergenceHistory - Sets the array used to hold the convergence history.
2582 
2583    Logically Collective on Tao
2584 
2585    Input Parameters:
2586 +  tao - the Tao solver context
2587 .  obj   - array to hold objective value history
2588 .  resid - array to hold residual history
2589 .  cnorm - array to hold constraint violation history
2590 .  lits - integer array holds the number of linear iterations for each Tao iteration
2591 .  na  - size of obj, resid, and cnorm
2592 -  reset - PetscTrue indicates each new minimization resets the history counter to zero,
2593            else it continues storing new values for new minimizations after the old ones
2594 
2595    Notes:
2596    If set, TAO will fill the given arrays with the indicated
2597    information at each iteration.  If 'obj','resid','cnorm','lits' are
2598    *all* NULL then space (using size na, or 1000 if na is PETSC_DECIDE or
2599    PETSC_DEFAULT) is allocated for the history.
2600    If not all are NULL, then only the non-NULL information categories
2601    will be stored, the others will be ignored.
2602 
2603    Any convergence information after iteration number 'na' will not be stored.
2604 
2605    This routine is useful, e.g., when running a code for purposes
2606    of accurate performance monitoring, when no I/O should be done
2607    during the section of code that is being timed.
2608 
2609    Level: intermediate
2610 
2611 .seealso: `TaoGetConvergenceHistory()`
2612 
2613 @*/
2614 PetscErrorCode TaoSetConvergenceHistory(Tao tao, PetscReal obj[], PetscReal resid[], PetscReal cnorm[], PetscInt lits[], PetscInt na,PetscBool reset)
2615 {
2616   PetscFunctionBegin;
2617   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2618   if (obj) PetscValidRealPointer(obj,2);
2619   if (resid) PetscValidRealPointer(resid,3);
2620   if (cnorm) PetscValidRealPointer(cnorm,4);
2621   if (lits) PetscValidIntPointer(lits,5);
2622 
2623   if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
2624   if (!obj && !resid && !cnorm && !lits) {
2625     PetscCall(PetscCalloc4(na,&obj,na,&resid,na,&cnorm,na,&lits));
2626     tao->hist_malloc = PETSC_TRUE;
2627   }
2628 
2629   tao->hist_obj = obj;
2630   tao->hist_resid = resid;
2631   tao->hist_cnorm = cnorm;
2632   tao->hist_lits = lits;
2633   tao->hist_max   = na;
2634   tao->hist_reset = reset;
2635   tao->hist_len = 0;
2636   PetscFunctionReturn(0);
2637 }
2638 
2639 /*@C
2640    TaoGetConvergenceHistory - Gets the arrays used to hold the convergence history.
2641 
2642    Collective on Tao
2643 
2644    Input Parameter:
2645 .  tao - the Tao context
2646 
2647    Output Parameters:
2648 +  obj   - array used to hold objective value history
2649 .  resid - array used to hold residual history
2650 .  cnorm - array used to hold constraint violation history
2651 .  lits  - integer array used to hold linear solver iteration count
2652 -  nhist  - size of obj, resid, cnorm, and lits
2653 
2654    Notes:
2655     This routine must be preceded by calls to TaoSetConvergenceHistory()
2656     and TaoSolve(), otherwise it returns useless information.
2657 
2658     The calling sequence for this routine in Fortran is
2659 $   call TaoGetConvergenceHistory(Tao tao, PetscInt nhist, PetscErrorCode ierr)
2660 
2661    This routine is useful, e.g., when running a code for purposes
2662    of accurate performance monitoring, when no I/O should be done
2663    during the section of code that is being timed.
2664 
2665    Level: advanced
2666 
2667 .seealso: `TaoSetConvergenceHistory()`
2668 
2669 @*/
2670 PetscErrorCode TaoGetConvergenceHistory(Tao tao, PetscReal **obj, PetscReal **resid, PetscReal **cnorm, PetscInt **lits, PetscInt *nhist)
2671 {
2672   PetscFunctionBegin;
2673   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2674   if (obj)   *obj   = tao->hist_obj;
2675   if (cnorm) *cnorm = tao->hist_cnorm;
2676   if (resid) *resid = tao->hist_resid;
2677   if (nhist) *nhist = tao->hist_len;
2678   PetscFunctionReturn(0);
2679 }
2680 
2681 /*@
2682    TaoSetApplicationContext - Sets the optional user-defined context for
2683    a solver.
2684 
2685    Logically Collective on Tao
2686 
2687    Input Parameters:
2688 +  tao  - the Tao context
2689 -  usrP - optional user context
2690 
2691    Level: intermediate
2692 
2693 .seealso: `TaoGetApplicationContext()`, `TaoSetApplicationContext()`
2694 @*/
2695 PetscErrorCode  TaoSetApplicationContext(Tao tao,void *usrP)
2696 {
2697   PetscFunctionBegin;
2698   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2699   tao->user = usrP;
2700   PetscFunctionReturn(0);
2701 }
2702 
2703 /*@
2704    TaoGetApplicationContext - Gets the user-defined context for a
2705    TAO solvers.
2706 
2707    Not Collective
2708 
2709    Input Parameter:
2710 .  tao  - Tao context
2711 
2712    Output Parameter:
2713 .  usrP - user context
2714 
2715    Level: intermediate
2716 
2717 .seealso: `TaoSetApplicationContext()`
2718 @*/
2719 PetscErrorCode  TaoGetApplicationContext(Tao tao,void *usrP)
2720 {
2721   PetscFunctionBegin;
2722   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2723   PetscValidPointer(usrP,2);
2724   *(void**)usrP = tao->user;
2725   PetscFunctionReturn(0);
2726 }
2727 
2728 /*@
2729    TaoSetGradientNorm - Sets the matrix used to define the inner product that measures the size of the gradient.
2730 
2731    Collective on tao
2732 
2733    Input Parameters:
2734 +  tao  - the Tao context
2735 -  M    - gradient norm
2736 
2737    Level: beginner
2738 
2739 .seealso: `TaoGetGradientNorm()`, `TaoGradientNorm()`
2740 @*/
2741 PetscErrorCode  TaoSetGradientNorm(Tao tao, Mat M)
2742 {
2743   PetscFunctionBegin;
2744   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2745   PetscValidHeaderSpecific(M,MAT_CLASSID,2);
2746   PetscCall(PetscObjectReference((PetscObject)M));
2747   PetscCall(MatDestroy(&tao->gradient_norm));
2748   PetscCall(VecDestroy(&tao->gradient_norm_tmp));
2749   tao->gradient_norm = M;
2750   PetscCall(MatCreateVecs(M, NULL, &tao->gradient_norm_tmp));
2751   PetscFunctionReturn(0);
2752 }
2753 
2754 /*@
2755    TaoGetGradientNorm - Returns the matrix used to define the inner product for measuring the size of the gradient.
2756 
2757    Not Collective
2758 
2759    Input Parameter:
2760 .  tao  - Tao context
2761 
2762    Output Parameter:
2763 .  M - gradient norm
2764 
2765    Level: beginner
2766 
2767 .seealso: `TaoSetGradientNorm()`, `TaoGradientNorm()`
2768 @*/
2769 PetscErrorCode  TaoGetGradientNorm(Tao tao, Mat *M)
2770 {
2771   PetscFunctionBegin;
2772   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2773   PetscValidPointer(M,2);
2774   *M = tao->gradient_norm;
2775   PetscFunctionReturn(0);
2776 }
2777 
2778 /*@C
2779    TaoGradientNorm - Compute the norm with respect to the inner product the user has set.
2780 
2781    Collective on tao
2782 
2783    Input Parameters:
2784 +  tao      - the Tao context
2785 .  gradient - the gradient to be computed
2786 -  norm     - the norm type
2787 
2788    Output Parameter:
2789 .  gnorm    - the gradient norm
2790 
2791    Level: developer
2792 
2793 .seealso: `TaoSetGradientNorm()`, `TaoGetGradientNorm()`
2794 @*/
2795 PetscErrorCode  TaoGradientNorm(Tao tao, Vec gradient, NormType type, PetscReal *gnorm)
2796 {
2797   PetscFunctionBegin;
2798   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2799   PetscValidHeaderSpecific(gradient,VEC_CLASSID,2);
2800   PetscValidLogicalCollectiveEnum(tao,type,3);
2801   PetscValidRealPointer(gnorm,4);
2802   if (tao->gradient_norm) {
2803     PetscScalar gnorms;
2804 
2805     PetscCheck(type == NORM_2,PetscObjectComm((PetscObject)gradient), PETSC_ERR_ARG_WRONG, "Norm type must be NORM_2 if an inner product for the gradient norm is set.");
2806     PetscCall(MatMult(tao->gradient_norm, gradient, tao->gradient_norm_tmp));
2807     PetscCall(VecDot(gradient, tao->gradient_norm_tmp, &gnorms));
2808     *gnorm = PetscRealPart(PetscSqrtScalar(gnorms));
2809   } else {
2810     PetscCall(VecNorm(gradient, type, gnorm));
2811   }
2812   PetscFunctionReturn(0);
2813 }
2814 
2815 /*@C
2816    TaoMonitorDrawCtxCreate - Creates the monitor context for TaoMonitorDrawCtx
2817 
2818    Collective on Tao
2819 
2820    Output Patameter:
2821 .    ctx - the monitor context
2822 
2823    Options Database:
2824 .   -tao_draw_solution_initial - show initial guess as well as current solution
2825 
2826    Level: intermediate
2827 
2828 .seealso: `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorDrawCtx()`
2829 @*/
2830 PetscErrorCode  TaoMonitorDrawCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TaoMonitorDrawCtx *ctx)
2831 {
2832   PetscFunctionBegin;
2833   PetscCall(PetscNew(ctx));
2834   PetscCall(PetscViewerDrawOpen(comm,host,label,x,y,m,n,&(*ctx)->viewer));
2835   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
2836   (*ctx)->howoften = howoften;
2837   PetscFunctionReturn(0);
2838 }
2839 
2840 /*@C
2841    TaoMonitorDrawCtxDestroy - Destroys the monitor context for TaoMonitorDrawSolution()
2842 
2843    Collective on Tao
2844 
2845    Input Parameters:
2846 .    ctx - the monitor context
2847 
2848    Level: intermediate
2849 
2850 .seealso: `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorDrawSolution()`
2851 @*/
2852 PetscErrorCode  TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx *ictx)
2853 {
2854   PetscFunctionBegin;
2855   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
2856   PetscCall(PetscFree(*ictx));
2857   PetscFunctionReturn(0);
2858 }
2859