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