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