xref: /petsc/src/tao/interface/taosolver.c (revision df4cd43f92eaa320656440c40edb1046daee8f75)
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, 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: [](chapter_tao), `Tao`, `TaoSolve()`, `TaoDestroy()`, `TAOSetFromOptions()`, `TAOSetType()`
87 @*/
88 PetscErrorCode TaoCreate(MPI_Comm comm, Tao *newtao)
89 {
90   Tao tao;
91 
92   PetscFunctionBegin;
93   PetscValidPointer(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: [](chapter_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: [](chapter_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: [](chapter_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    Reference:
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: [](chapter_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: [](chapter_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: [](chapter_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   InputParameters:
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: [](chapter_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: [](chapter_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: [](chapter_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   PetscValidBoolPointer(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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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    OutputParameter:
957 .  fmin - the minimum function value
958 
959    Level: intermediate
960 
961 .seealso: [](chapter_tao), `Tao`, `TaoConvergedReason`, `TaoSetFunctionLowerBound()`
962 @*/
963 PetscErrorCode TaoGetFunctionLowerBound(Tao tao, PetscReal *fmin)
964 {
965   PetscFunctionBegin;
966   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
967   PetscValidRealPointer(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: [](chapter_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: [](chapter_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1016 @*/
1017 PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao tao, PetscInt *nfcn)
1018 {
1019   PetscFunctionBegin;
1020   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1021   PetscValidIntPointer(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: [](chapter_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1040 @*/
1041 PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao tao, PetscInt *nfuncs)
1042 {
1043   PetscFunctionBegin;
1044   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1045   PetscValidIntPointer(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: [](chapter_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: [](chapter_tao), `Tao`, `TaoSetMaximumIterations()`, `TaoGetMaximumFunctionEvaluations()`
1090 @*/
1091 PetscErrorCode TaoGetMaximumIterations(Tao tao, PetscInt *maxits)
1092 {
1093   PetscFunctionBegin;
1094   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1095   PetscValidIntPointer(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: [](chapter_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: [](chapter_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetCurrentTrustRegionRadius()`, `TAONTR`
1140 @*/
1141 PetscErrorCode TaoGetInitialTrustRegionRadius(Tao tao, PetscReal *radius)
1142 {
1143   PetscFunctionBegin;
1144   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1145   PetscValidRealPointer(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: [](chapter_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetInitialTrustRegionRadius()`, `TAONTR`
1164 @*/
1165 PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao tao, PetscReal *radius)
1166 {
1167   PetscFunctionBegin;
1168   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1169   PetscValidRealPointer(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: [](chapter_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: [](chapter_tao), `Tao`, `KSP`
1218 @*/
1219 PetscErrorCode TaoGetKSP(Tao tao, KSP *ksp)
1220 {
1221   PetscFunctionBegin;
1222   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1223   PetscValidPointer(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: [](chapter_tao), `Tao`, `TaoGetKSP()`
1246 @*/
1247 PetscErrorCode TaoGetLinearSolveIterations(Tao tao, PetscInt *lits)
1248 {
1249   PetscFunctionBegin;
1250   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1251   PetscValidIntPointer(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: [](chapter_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`
1270 @*/
1271 PetscErrorCode TaoGetLineSearch(Tao tao, TaoLineSearch *ls)
1272 {
1273   PetscFunctionBegin;
1274   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1275   PetscValidPointer(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 -  ls - the line search used in the optimization solver
1287 
1288    Level: developer
1289 
1290 .seealso: [](chapter_tao), `Tao`, `TaoGetLineSearch()`, `TaoLineSearchApply()`
1291 @*/
1292 PetscErrorCode TaoAddLineSearchCounts(Tao tao)
1293 {
1294   PetscBool flg;
1295   PetscInt  nfeval, ngeval, nfgeval;
1296 
1297   PetscFunctionBegin;
1298   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1299   if (tao->linesearch) {
1300     PetscCall(TaoLineSearchIsUsingTaoRoutines(tao->linesearch, &flg));
1301     if (!flg) {
1302       PetscCall(TaoLineSearchGetNumberFunctionEvaluations(tao->linesearch, &nfeval, &ngeval, &nfgeval));
1303       tao->nfuncs += nfeval;
1304       tao->ngrads += ngeval;
1305       tao->nfuncgrads += nfgeval;
1306     }
1307   }
1308   PetscFunctionReturn(PETSC_SUCCESS);
1309 }
1310 
1311 /*@
1312   TaoGetSolution - Returns the vector with the current solution from the `Tao` object
1313 
1314   Not Collective
1315 
1316   Input Parameter:
1317 . tao - the `Tao` context
1318 
1319   Output Parameter:
1320 . X - the current solution
1321 
1322   Level: intermediate
1323 
1324   Note:
1325   The returned vector will be the same object that was passed into `TaoSetSolution()`
1326 
1327 .seealso: [](chapter_tao), `Tao`, `TaoSetSolution()`, `TaoSolve()`
1328 @*/
1329 PetscErrorCode TaoGetSolution(Tao tao, Vec *X)
1330 {
1331   PetscFunctionBegin;
1332   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1333   PetscValidPointer(X, 2);
1334   *X = tao->solution;
1335   PetscFunctionReturn(PETSC_SUCCESS);
1336 }
1337 
1338 /*@
1339    TaoResetStatistics - Initialize the statistics collected by the `Tao` object.
1340    These statistics include the iteration number, residual norms, and convergence status.
1341    This routine gets called before solving each optimization problem.
1342 
1343    Collective
1344 
1345    Input Parameter:
1346 .  solver - the `Tao` context
1347 
1348    Level: developer
1349 
1350 .seealso: [](chapter_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
1351 @*/
1352 PetscErrorCode TaoResetStatistics(Tao tao)
1353 {
1354   PetscFunctionBegin;
1355   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1356   tao->niter        = 0;
1357   tao->nfuncs       = 0;
1358   tao->nfuncgrads   = 0;
1359   tao->ngrads       = 0;
1360   tao->nhess        = 0;
1361   tao->njac         = 0;
1362   tao->nconstraints = 0;
1363   tao->ksp_its      = 0;
1364   tao->ksp_tot_its  = 0;
1365   tao->reason       = TAO_CONTINUE_ITERATING;
1366   tao->residual     = 0.0;
1367   tao->cnorm        = 0.0;
1368   tao->step         = 0.0;
1369   tao->lsflag       = PETSC_FALSE;
1370   if (tao->hist_reset) tao->hist_len = 0;
1371   PetscFunctionReturn(PETSC_SUCCESS);
1372 }
1373 
1374 /*@C
1375   TaoSetUpdate - Sets the general-purpose update function called
1376   at the beginning of every iteration of the optimization algorithm. Specifically
1377   it is called at the top of every iteration, after the new solution and the gradient
1378   is determined, but before the Hessian is computed (if applicable).
1379 
1380   Logically Collective
1381 
1382   Input Parameters:
1383 + tao - The tao solver context
1384 - func - The function
1385 
1386   Calling sequence of func:
1387 $ func (Tao tao, PetscInt step);
1388 
1389 . step - The current step of the iteration
1390 
1391   Level: advanced
1392 
1393 .seealso: [](chapter_tao), `Tao`, `TaoSolve()`
1394 @*/
1395 PetscErrorCode TaoSetUpdate(Tao tao, PetscErrorCode (*func)(Tao, PetscInt, void *), void *ctx)
1396 {
1397   PetscFunctionBegin;
1398   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1399   tao->ops->update = func;
1400   tao->user_update = ctx;
1401   PetscFunctionReturn(PETSC_SUCCESS);
1402 }
1403 
1404 /*@C
1405   TaoSetConvergenceTest - Sets the function that is to be used to test
1406   for convergence o fthe iterative minimization solution.  The new convergence
1407   testing routine will replace Tao's default convergence test.
1408 
1409   Logically Collective
1410 
1411   Input Parameters:
1412 + tao - the `Tao` object
1413 . conv - the routine to test for convergence
1414 - ctx - [optional] context for private data for the convergence routine
1415         (may be `NULL`)
1416 
1417   Calling sequence of conv:
1418 $   PetscErrorCode conv(Tao tao, void *ctx)
1419 
1420 + tao - the `Tao` object
1421 - ctx - [optional] convergence context
1422 
1423   Level: advanced
1424 
1425   Note:
1426   The new convergence testing routine should call `TaoSetConvergedReason()`.
1427 
1428 .seealso: [](chapter_tao), `Tao`, `TaoSolve()`, `TaoSetConvergedReason()`, `TaoGetSolutionStatus()`, `TaoGetTolerances()`, `TaoSetMonitor`
1429 @*/
1430 PetscErrorCode TaoSetConvergenceTest(Tao tao, PetscErrorCode (*conv)(Tao, void *), void *ctx)
1431 {
1432   PetscFunctionBegin;
1433   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1434   tao->ops->convergencetest = conv;
1435   tao->cnvP                 = ctx;
1436   PetscFunctionReturn(PETSC_SUCCESS);
1437 }
1438 
1439 /*@C
1440    TaoSetMonitor - Sets an additional function that is to be used at every
1441    iteration of the solver to display the iteration's
1442    progress.
1443 
1444    Logically Collective
1445 
1446    Input Parameters:
1447 +  tao - the `Tao` solver context
1448 .  mymonitor - monitoring routine
1449 -  mctx - [optional] user-defined context for private data for the
1450           monitor routine (may be `NULL`)
1451 
1452    Calling sequence of mymonitor:
1453 .vb
1454      PetscErrorCode mymonitor(Tao tao,void *mctx)
1455 .ve
1456 
1457 +    tao - the `Tao` solver context
1458 -    mctx - [optional] monitoring context
1459 
1460    Options Database Keys:
1461 +    -tao_monitor        - sets the default monitor `TaoMonitorDefault()`
1462 .    -tao_smonitor       - sets short monitor
1463 .    -tao_cmonitor       - same as smonitor plus constraint norm
1464 .    -tao_view_solution   - view solution at each iteration
1465 .    -tao_view_gradient   - view gradient at each iteration
1466 .    -tao_view_ls_residual - view least-squares residual vector at each iteration
1467 -    -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.
1468 
1469    Level: intermediate
1470 
1471    Notes:
1472    Several different monitoring routines may be set by calling
1473    `TaoSetMonitor()` multiple times; all will be called in the
1474    order in which they were set.
1475 
1476    Fortran Note:
1477     Only one monitor function may be set
1478 
1479 .seealso: [](chapter_tao), `Tao`, `TaoSolve()`, `TaoMonitorDefault()`, `TaoCancelMonitors()`, `TaoSetDestroyRoutine()`, `TaoView()`
1480 @*/
1481 PetscErrorCode TaoSetMonitor(Tao tao, PetscErrorCode (*func)(Tao, void *), void *ctx, PetscErrorCode (*dest)(void **))
1482 {
1483   PetscInt  i;
1484   PetscBool identical;
1485 
1486   PetscFunctionBegin;
1487   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1488   PetscCheck(tao->numbermonitors < MAXTAOMONITORS, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Cannot attach another monitor -- max=%d", MAXTAOMONITORS);
1489 
1490   for (i = 0; i < tao->numbermonitors; i++) {
1491     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))func, ctx, dest, (PetscErrorCode(*)(void))tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i], &identical));
1492     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
1493   }
1494   tao->monitor[tao->numbermonitors]        = func;
1495   tao->monitorcontext[tao->numbermonitors] = (void *)ctx;
1496   tao->monitordestroy[tao->numbermonitors] = dest;
1497   ++tao->numbermonitors;
1498   PetscFunctionReturn(PETSC_SUCCESS);
1499 }
1500 
1501 /*@
1502    TaoCancelMonitors - Clears all the monitor functions for a `Tao` object.
1503 
1504    Logically Collective
1505 
1506    Input Parameter:
1507 .  tao - the `Tao` solver context
1508 
1509    Options Database Key:
1510 .  -tao_cancelmonitors - cancels all monitors that have been hardwired
1511     into a code by calls to `TaoSetMonitor()`, but does not cancel those
1512     set via the options database
1513 
1514    Level: advanced
1515 
1516    Note:
1517    There is no way to clear one specific monitor from a `Tao` object.
1518 
1519 .seealso: [](chapter_tao), `Tao`, `TaoMonitorDefault()`, `TaoSetMonitor()`
1520 @*/
1521 PetscErrorCode TaoCancelMonitors(Tao tao)
1522 {
1523   PetscInt i;
1524 
1525   PetscFunctionBegin;
1526   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1527   for (i = 0; i < tao->numbermonitors; i++) {
1528     if (tao->monitordestroy[i]) PetscCall((*tao->monitordestroy[i])(&tao->monitorcontext[i]));
1529   }
1530   tao->numbermonitors = 0;
1531   PetscFunctionReturn(PETSC_SUCCESS);
1532 }
1533 
1534 /*@
1535    TaoMonitorDefault - Default routine for monitoring progress of `TaoSolve()`
1536 
1537    Collective
1538 
1539    Input Parameters:
1540 +  tao - the `Tao` context
1541 -  ctx - `PetscViewer` context or `NULL`
1542 
1543    Options Database Key:
1544 .  -tao_monitor - turn on default monitoring
1545 
1546    Level: advanced
1547 
1548    Note:
1549    This monitor prints the function value and gradient
1550    norm at each iteration.
1551 
1552 .seealso: [](chapter_tao), `Tao`, `TaoDefaultSMonitor()`, `TaoSetMonitor()`
1553 @*/
1554 PetscErrorCode TaoMonitorDefault(Tao tao, void *ctx)
1555 {
1556   PetscInt    its, tabs;
1557   PetscReal   fct, gnorm;
1558   PetscViewer viewer = (PetscViewer)ctx;
1559 
1560   PetscFunctionBegin;
1561   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1562   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1563   its   = tao->niter;
1564   fct   = tao->fc;
1565   gnorm = tao->residual;
1566   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1567   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1568   if (its == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
1569     PetscCall(PetscViewerASCIIPrintf(viewer, "  Iteration information for %s solve.\n", ((PetscObject)tao)->prefix));
1570     tao->header_printed = PETSC_TRUE;
1571   }
1572   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " TAO,", its));
1573   PetscCall(PetscViewerASCIIPrintf(viewer, "  Function value: %g,", (double)fct));
1574   if (gnorm >= PETSC_INFINITY) {
1575     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: Inf \n"));
1576   } else {
1577     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g \n", (double)gnorm));
1578   }
1579   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1580   PetscFunctionReturn(PETSC_SUCCESS);
1581 }
1582 
1583 /*@
1584    TaoDefaultGMonitor - Default routine for monitoring progress of `TaoSolve()` with extra detail on the globalization method.
1585 
1586    Collective
1587 
1588    Input Parameters:
1589 +  tao - the `Tao` context
1590 -  ctx - `PetscViewer` context or `NULL`
1591 
1592    Options Database Key:
1593 .  -tao_gmonitor - turn on monitoring with globalization information
1594 
1595    Level: advanced
1596 
1597    Note:
1598    This monitor prints the function value and gradient norm at each
1599    iteration, as well as the step size and trust radius. Note that the
1600    step size and trust radius may be the same for some algorithms.
1601 
1602 .seealso: [](chapter_tao), `Tao`, `TaoDefaultSMonitor()`, `TaoSetMonitor()`
1603 @*/
1604 PetscErrorCode TaoDefaultGMonitor(Tao tao, void *ctx)
1605 {
1606   PetscInt    its, tabs;
1607   PetscReal   fct, gnorm, stp, tr;
1608   PetscViewer viewer = (PetscViewer)ctx;
1609 
1610   PetscFunctionBegin;
1611   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1612   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1613   its   = tao->niter;
1614   fct   = tao->fc;
1615   gnorm = tao->residual;
1616   stp   = tao->step;
1617   tr    = tao->trust;
1618   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1619   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1620   if (its == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
1621     PetscCall(PetscViewerASCIIPrintf(viewer, "  Iteration information for %s solve.\n", ((PetscObject)tao)->prefix));
1622     tao->header_printed = PETSC_TRUE;
1623   }
1624   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " TAO,", its));
1625   PetscCall(PetscViewerASCIIPrintf(viewer, "  Function value: %g,", (double)fct));
1626   if (gnorm >= PETSC_INFINITY) {
1627     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: Inf,"));
1628   } else {
1629     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g,", (double)gnorm));
1630   }
1631   PetscCall(PetscViewerASCIIPrintf(viewer, "  Step: %g,  Trust: %g\n", (double)stp, (double)tr));
1632   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1633   PetscFunctionReturn(PETSC_SUCCESS);
1634 }
1635 
1636 /*@
1637    TaoDefaultSMonitor - Default routine for monitoring progress of `TaoSolve()`
1638 
1639    Collective
1640 
1641    Input Parameters:
1642 +  tao - the `Tao` context
1643 -  ctx - `PetscViewer` context of type `PETSCVIEWERASCII`
1644 
1645    Options Database Key:
1646 .  -tao_smonitor - turn on default short monitoring
1647 
1648    Level: advanced
1649 
1650    Note:
1651    Same as `TaoMonitorDefault()` except
1652    it prints fewer digits of the residual as the residual gets smaller.
1653    This is because the later digits are meaningless and are often
1654    different on different machines; by using this routine different
1655    machines will usually generate the same output.
1656 
1657 .seealso: [](chapter_tao), `Tao`, `TaoMonitorDefault()`, `TaoSetMonitor()`
1658 @*/
1659 PetscErrorCode TaoDefaultSMonitor(Tao tao, void *ctx)
1660 {
1661   PetscInt    its, tabs;
1662   PetscReal   fct, gnorm;
1663   PetscViewer viewer = (PetscViewer)ctx;
1664 
1665   PetscFunctionBegin;
1666   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1667   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1668   its   = tao->niter;
1669   fct   = tao->fc;
1670   gnorm = tao->residual;
1671   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1672   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1673   PetscCall(PetscViewerASCIIPrintf(viewer, "iter = %3" PetscInt_FMT ",", its));
1674   PetscCall(PetscViewerASCIIPrintf(viewer, " Function value %g,", (double)fct));
1675   if (gnorm >= PETSC_INFINITY) {
1676     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: Inf \n"));
1677   } else if (gnorm > 1.e-6) {
1678     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: %g \n", (double)gnorm));
1679   } else if (gnorm > 1.e-11) {
1680     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: < 1.0e-6 \n"));
1681   } else {
1682     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: < 1.0e-11 \n"));
1683   }
1684   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1685   PetscFunctionReturn(PETSC_SUCCESS);
1686 }
1687 
1688 /*@
1689    TaoDefaultCMonitor - same as `TaoMonitorDefault()` except
1690    it prints the norm of the constraint function.
1691 
1692    Collective
1693 
1694    Input Parameters:
1695 +  tao - the `Tao` context
1696 -  ctx - `PetscViewer` context or `NULL`
1697 
1698    Options Database Key:
1699 .  -tao_cmonitor - monitor the constraints
1700 
1701    Level: advanced
1702 
1703 .seealso: [](chapter_tao), `Tao`, `TaoMonitorDefault()`, `TaoSetMonitor()`
1704 @*/
1705 PetscErrorCode TaoDefaultCMonitor(Tao tao, void *ctx)
1706 {
1707   PetscInt    its, tabs;
1708   PetscReal   fct, gnorm;
1709   PetscViewer viewer = (PetscViewer)ctx;
1710 
1711   PetscFunctionBegin;
1712   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1713   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1714   its   = tao->niter;
1715   fct   = tao->fc;
1716   gnorm = tao->residual;
1717   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1718   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1719   PetscCall(PetscViewerASCIIPrintf(viewer, "iter = %" PetscInt_FMT ",", its));
1720   PetscCall(PetscViewerASCIIPrintf(viewer, " Function value: %g,", (double)fct));
1721   PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g ", (double)gnorm));
1722   PetscCall(PetscViewerASCIIPrintf(viewer, "  Constraint: %g \n", (double)tao->cnorm));
1723   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1724   PetscFunctionReturn(PETSC_SUCCESS);
1725 }
1726 
1727 /*@C
1728    TaoSolutionMonitor - Views the solution at each iteration of `TaoSolve()`
1729 
1730    Collective
1731 
1732    Input Parameters:
1733 +  tao - the `Tao` context
1734 -  ctx - `PetscViewer` context or `NULL`
1735 
1736    Options Database Key:
1737 .  -tao_view_solution - view the solution
1738 
1739    Level: advanced
1740 
1741 .seealso: [](chapter_tao), `Tao`, `TaoDefaultSMonitor()`, `TaoSetMonitor()`
1742 @*/
1743 PetscErrorCode TaoSolutionMonitor(Tao tao, void *ctx)
1744 {
1745   PetscViewer viewer = (PetscViewer)ctx;
1746 
1747   PetscFunctionBegin;
1748   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1749   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1750   PetscCall(VecView(tao->solution, viewer));
1751   PetscFunctionReturn(PETSC_SUCCESS);
1752 }
1753 
1754 /*@C
1755    TaoGradientMonitor - Views the gradient at each iteration of `TaoSolve()`
1756 
1757    Collective
1758 
1759    Input Parameters:
1760 +  tao - the `Tao` context
1761 -  ctx - `PetscViewer` context or `NULL`
1762 
1763    Options Database Key:
1764 .  -tao_view_gradient - view the gradient at each iteration
1765 
1766    Level: advanced
1767 
1768 .seealso: [](chapter_tao), `Tao`, `TaoDefaultSMonitor()`, `TaoSetMonitor()`
1769 @*/
1770 PetscErrorCode TaoGradientMonitor(Tao tao, void *ctx)
1771 {
1772   PetscViewer viewer = (PetscViewer)ctx;
1773 
1774   PetscFunctionBegin;
1775   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1776   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1777   PetscCall(VecView(tao->gradient, viewer));
1778   PetscFunctionReturn(PETSC_SUCCESS);
1779 }
1780 
1781 /*@C
1782    TaoStepDirectionMonitor - Views the step-direction at each iteration of `TaoSolve()`
1783 
1784    Collective
1785 
1786    Input Parameters:
1787 +  tao - the `Tao` context
1788 -  ctx - `PetscViewer` context or `NULL`
1789 
1790    Options Database Key:
1791 .  -tao_view_stepdirection - view the step direction vector at each iteration
1792 
1793    Level: advanced
1794 
1795 .seealso: [](chapter_tao), `Tao`, `TaoDefaultSMonitor()`, `TaoSetMonitor()`
1796 @*/
1797 PetscErrorCode TaoStepDirectionMonitor(Tao tao, void *ctx)
1798 {
1799   PetscViewer viewer = (PetscViewer)ctx;
1800 
1801   PetscFunctionBegin;
1802   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1803   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1804   PetscCall(VecView(tao->stepdirection, viewer));
1805   PetscFunctionReturn(PETSC_SUCCESS);
1806 }
1807 
1808 /*@C
1809    TaoDrawSolutionMonitor - Plots the solution at each iteration of `TaoSolve()`
1810 
1811    Collective
1812 
1813    Input Parameters:
1814 +  tao - the `Tao` context
1815 -  ctx - `TaoMonitorDraw` context
1816 
1817    Options Database Key:
1818 .  -tao_draw_solution - draw the solution at each iteration
1819 
1820    Level: advanced
1821 
1822 .seealso: [](chapter_tao), `Tao`, `TaoSolutionMonitor()`, `TaoSetMonitor()`, `TaoDrawGradientMonitor`, `TaoMonitorDraw`
1823 @*/
1824 PetscErrorCode TaoDrawSolutionMonitor(Tao tao, void *ctx)
1825 {
1826   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;
1827 
1828   PetscFunctionBegin;
1829   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1830   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(PETSC_SUCCESS);
1831   PetscCall(VecView(tao->solution, ictx->viewer));
1832   PetscFunctionReturn(PETSC_SUCCESS);
1833 }
1834 
1835 /*@C
1836    TaoDrawGradientMonitor - Plots the gradient at each iteration of `TaoSolve()`
1837 
1838    Collective
1839 
1840    Input Parameters:
1841 +  tao - the `Tao` context
1842 -  ctx - `PetscViewer` context
1843 
1844    Options Database Key:
1845 .  -tao_draw_gradient - draw the gradient at each iteration
1846 
1847    Level: advanced
1848 
1849 .seealso: [](chapter_tao), `Tao`, `TaoGradientMonitor()`, `TaoSetMonitor()`, `TaoDrawSolutionMonitor`
1850 @*/
1851 PetscErrorCode TaoDrawGradientMonitor(Tao tao, void *ctx)
1852 {
1853   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;
1854 
1855   PetscFunctionBegin;
1856   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1857   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(PETSC_SUCCESS);
1858   PetscCall(VecView(tao->gradient, ictx->viewer));
1859   PetscFunctionReturn(PETSC_SUCCESS);
1860 }
1861 
1862 /*@C
1863    TaoDrawStepMonitor - Plots the step direction at each iteration of `TaoSolve()`
1864 
1865    Collective
1866 
1867    Input Parameters:
1868 +  tao - the `Tao` context
1869 -  ctx - the `PetscViewer` context
1870 
1871    Options Database Key:
1872 .  -tao_draw_step - draw the step direction at each iteration
1873 
1874    Level: advanced
1875 
1876 .seealso: [](chapter_tao), `Tao`, `TaoSetMonitor()`, `TaoDrawSolutionMonitor`
1877 @*/
1878 PetscErrorCode TaoDrawStepMonitor(Tao tao, void *ctx)
1879 {
1880   PetscViewer viewer = (PetscViewer)ctx;
1881 
1882   PetscFunctionBegin;
1883   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1884   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1885   PetscCall(VecView(tao->stepdirection, viewer));
1886   PetscFunctionReturn(PETSC_SUCCESS);
1887 }
1888 
1889 /*@C
1890    TaoResidualMonitor - Views the least-squares residual at each iteration of `TaoSolve()`
1891 
1892    Collective
1893 
1894    Input Parameters:
1895 +  tao - the `Tao` context
1896 -  ctx - the `PetscViewer` context or `NULL`
1897 
1898    Options Database Key:
1899 .  -tao_view_ls_residual - view the residual at each iteration
1900 
1901    Level: advanced
1902 
1903 .seealso: [](chapter_tao), `Tao`, `TaoDefaultSMonitor()`, `TaoSetMonitor()`
1904 @*/
1905 PetscErrorCode TaoResidualMonitor(Tao tao, void *ctx)
1906 {
1907   PetscViewer viewer = (PetscViewer)ctx;
1908 
1909   PetscFunctionBegin;
1910   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1911   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1912   PetscCall(VecView(tao->ls_res, viewer));
1913   PetscFunctionReturn(PETSC_SUCCESS);
1914 }
1915 
1916 /*@
1917    TaoDefaultConvergenceTest - Determines whether the solver should continue iterating
1918    or terminate.
1919 
1920    Collective
1921 
1922    Input Parameters:
1923 +  tao - the `Tao` context
1924 -  dummy - unused dummy context
1925 
1926    Output Parameter:
1927 .  reason - for terminating
1928 
1929    Level: developer
1930 
1931    Notes:
1932    This routine checks the residual in the optimality conditions, the
1933    relative residual in the optimity conditions, the number of function
1934    evaluations, and the function value to test convergence.  Some
1935    solvers may use different convergence routines.
1936 
1937 .seealso: [](chapter_tao), `Tao`, `TaoSetTolerances()`, `TaoGetConvergedReason()`, `TaoSetConvergedReason()`
1938 @*/
1939 PetscErrorCode TaoDefaultConvergenceTest(Tao tao, void *dummy)
1940 {
1941   PetscInt           niter = tao->niter, nfuncs = PetscMax(tao->nfuncs, tao->nfuncgrads);
1942   PetscInt           max_funcs = tao->max_funcs;
1943   PetscReal          gnorm = tao->residual, gnorm0 = tao->gnorm0;
1944   PetscReal          f = tao->fc, steptol = tao->steptol, trradius = tao->step;
1945   PetscReal          gatol = tao->gatol, grtol = tao->grtol, gttol = tao->gttol;
1946   PetscReal          catol = tao->catol, crtol = tao->crtol;
1947   PetscReal          fmin = tao->fmin, cnorm = tao->cnorm;
1948   TaoConvergedReason reason = tao->reason;
1949 
1950   PetscFunctionBegin;
1951   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1952   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
1953 
1954   if (PetscIsInfOrNanReal(f)) {
1955     PetscCall(PetscInfo(tao, "Failed to converged, function value is Inf or NaN\n"));
1956     reason = TAO_DIVERGED_NAN;
1957   } else if (f <= fmin && cnorm <= catol) {
1958     PetscCall(PetscInfo(tao, "Converged due to function value %g < minimum function value %g\n", (double)f, (double)fmin));
1959     reason = TAO_CONVERGED_MINF;
1960   } else if (gnorm <= gatol && cnorm <= catol) {
1961     PetscCall(PetscInfo(tao, "Converged due to residual norm ||g(X)||=%g < %g\n", (double)gnorm, (double)gatol));
1962     reason = TAO_CONVERGED_GATOL;
1963   } else if (f != 0 && PetscAbsReal(gnorm / f) <= grtol && cnorm <= crtol) {
1964     PetscCall(PetscInfo(tao, "Converged due to residual ||g(X)||/|f(X)| =%g < %g\n", (double)(gnorm / f), (double)grtol));
1965     reason = TAO_CONVERGED_GRTOL;
1966   } else if (gnorm0 != 0 && ((gttol == 0 && gnorm == 0) || gnorm / gnorm0 < gttol) && cnorm <= crtol) {
1967     PetscCall(PetscInfo(tao, "Converged due to relative residual norm ||g(X)||/||g(X0)|| = %g < %g\n", (double)(gnorm / gnorm0), (double)gttol));
1968     reason = TAO_CONVERGED_GTTOL;
1969   } else if (max_funcs >= 0 && nfuncs > max_funcs) {
1970     PetscCall(PetscInfo(tao, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", nfuncs, max_funcs));
1971     reason = TAO_DIVERGED_MAXFCN;
1972   } else if (tao->lsflag != 0) {
1973     PetscCall(PetscInfo(tao, "Tao Line Search failure.\n"));
1974     reason = TAO_DIVERGED_LS_FAILURE;
1975   } else if (trradius < steptol && niter > 0) {
1976     PetscCall(PetscInfo(tao, "Trust region/step size too small: %g < %g\n", (double)trradius, (double)steptol));
1977     reason = TAO_CONVERGED_STEPTOL;
1978   } else if (niter >= tao->max_it) {
1979     PetscCall(PetscInfo(tao, "Exceeded maximum number of iterations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", niter, tao->max_it));
1980     reason = TAO_DIVERGED_MAXITS;
1981   } else {
1982     reason = TAO_CONTINUE_ITERATING;
1983   }
1984   tao->reason = reason;
1985   PetscFunctionReturn(PETSC_SUCCESS);
1986 }
1987 
1988 /*@C
1989    TaoSetOptionsPrefix - Sets the prefix used for searching for all
1990    Tao options in the database.
1991 
1992    Logically Collective
1993 
1994    Input Parameters:
1995 +  tao - the `Tao` context
1996 -  prefix - the prefix string to prepend to all Tao option requests
1997 
1998    Notes:
1999    A hyphen (-) must NOT be given at the beginning of the prefix name.
2000    The first character of all runtime options is AUTOMATICALLY the hyphen.
2001 
2002    For example, to distinguish between the runtime options for two
2003    different Tao solvers, one could call
2004 .vb
2005       TaoSetOptionsPrefix(tao1,"sys1_")
2006       TaoSetOptionsPrefix(tao2,"sys2_")
2007 .ve
2008 
2009    This would enable use of different options for each system, such as
2010 .vb
2011       -sys1_tao_method blmvm -sys1_tao_grtol 1.e-3
2012       -sys2_tao_method lmvm  -sys2_tao_grtol 1.e-4
2013 .ve
2014 
2015    Level: advanced
2016 
2017 .seealso: [](chapter_tao), `Tao`, `TaoSetFromOptions()`, `TaoAppendOptionsPrefix()`, `TaoGetOptionsPrefix()`
2018 @*/
2019 PetscErrorCode TaoSetOptionsPrefix(Tao tao, const char p[])
2020 {
2021   PetscFunctionBegin;
2022   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2023   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tao, p));
2024   if (tao->linesearch) PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, p));
2025   if (tao->ksp) PetscCall(KSPSetOptionsPrefix(tao->ksp, p));
2026   PetscFunctionReturn(PETSC_SUCCESS);
2027 }
2028 
2029 /*@C
2030    TaoAppendOptionsPrefix - Appends to the prefix used for searching for all Tao options in the database.
2031 
2032    Logically Collective
2033 
2034    Input Parameters:
2035 +  tao - the `Tao` solver context
2036 -  prefix - the prefix string to prepend to all `Tao` option requests
2037 
2038    Note:
2039    A hyphen (-) must NOT be given at the beginning of the prefix name.
2040    The first character of all runtime options is automatically the hyphen.
2041 
2042    Level: advanced
2043 
2044 .seealso: [](chapter_tao), `Tao`, `TaoSetFromOptions()`, `TaoSetOptionsPrefix()`, `TaoGetOptionsPrefix()`
2045 @*/
2046 PetscErrorCode TaoAppendOptionsPrefix(Tao tao, const char p[])
2047 {
2048   PetscFunctionBegin;
2049   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2050   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao, p));
2051   if (tao->linesearch) PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao->linesearch, p));
2052   if (tao->ksp) PetscCall(KSPAppendOptionsPrefix(tao->ksp, p));
2053   PetscFunctionReturn(PETSC_SUCCESS);
2054 }
2055 
2056 /*@C
2057   TaoGetOptionsPrefix - Gets the prefix used for searching for all
2058   Tao options in the database
2059 
2060   Not Collective
2061 
2062   Input Parameter:
2063 . tao - the `Tao` context
2064 
2065   Output Parameter:
2066 . prefix - pointer to the prefix string used is returned
2067 
2068   Fortran Note:
2069   Pass in a string 'prefix' of sufficient length to hold the prefix.
2070 
2071   Level: advanced
2072 
2073 .seealso: [](chapter_tao), `Tao`, `TaoSetFromOptions()`, `TaoSetOptionsPrefix()`, `TaoAppendOptionsPrefix()`
2074 @*/
2075 PetscErrorCode TaoGetOptionsPrefix(Tao tao, const char *p[])
2076 {
2077   PetscFunctionBegin;
2078   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2079   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)tao, p));
2080   PetscFunctionReturn(PETSC_SUCCESS);
2081 }
2082 
2083 /*@C
2084    TaoSetType - Sets the `TaoType` for the minimization solver.
2085 
2086    Collective
2087 
2088    Input Parameters:
2089 +  solver - the `Tao` solver context
2090 -  type - a known method
2091 
2092    Options Database Key:
2093 .  -tao_type <type> - Sets the method; use -help for a list
2094    of available methods (for instance, "-tao_type lmvm" or "-tao_type tron")
2095 
2096   Level: intermediate
2097 
2098 .seealso: [](chapter_tao), `Tao`, `TaoCreate()`, `TaoGetType()`, `TaoType`
2099 @*/
2100 PetscErrorCode TaoSetType(Tao tao, TaoType type)
2101 {
2102   PetscErrorCode (*create_xxx)(Tao);
2103   PetscBool issame;
2104 
2105   PetscFunctionBegin;
2106   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2107 
2108   PetscCall(PetscObjectTypeCompare((PetscObject)tao, type, &issame));
2109   if (issame) PetscFunctionReturn(PETSC_SUCCESS);
2110 
2111   PetscCall(PetscFunctionListFind(TaoList, type, (void (**)(void)) & create_xxx));
2112   PetscCheck(create_xxx, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Tao type %s", type);
2113 
2114   /* Destroy the existing solver information */
2115   PetscTryTypeMethod(tao, destroy);
2116   PetscCall(KSPDestroy(&tao->ksp));
2117   PetscCall(TaoLineSearchDestroy(&tao->linesearch));
2118   tao->ops->setup          = NULL;
2119   tao->ops->solve          = NULL;
2120   tao->ops->view           = NULL;
2121   tao->ops->setfromoptions = NULL;
2122   tao->ops->destroy        = NULL;
2123 
2124   tao->setupcalled = PETSC_FALSE;
2125 
2126   PetscCall((*create_xxx)(tao));
2127   PetscCall(PetscObjectChangeTypeName((PetscObject)tao, type));
2128   PetscFunctionReturn(PETSC_SUCCESS);
2129 }
2130 
2131 /*@C
2132    TaoRegister - Adds a method to the Tao package for minimization.
2133 
2134    Synopsis:
2135    TaoRegister(char *name_solver,char *path,char *name_Create,PetscErrorCode (*routine_Create)(Tao))
2136 
2137    Not collective
2138 
2139    Input Parameters:
2140 +  sname - name of a new user-defined solver
2141 -  func - routine to Create method context
2142 
2143    Sample usage:
2144 .vb
2145    TaoRegister("my_solver",MySolverCreate);
2146 .ve
2147 
2148    Then, your solver can be chosen with the procedural interface via
2149 $     TaoSetType(tao,"my_solver")
2150    or at runtime via the option
2151 $     -tao_type my_solver
2152 
2153    Level: advanced
2154 
2155    Note:
2156    `TaoRegister()` may be called multiple times to add several user-defined solvers.
2157 
2158 .seealso: [](chapter_tao), `Tao`, `TaoSetType()`, `TaoRegisterAll()`, `TaoRegisterDestroy()`
2159 @*/
2160 PetscErrorCode TaoRegister(const char sname[], PetscErrorCode (*func)(Tao))
2161 {
2162   PetscFunctionBegin;
2163   PetscCall(TaoInitializePackage());
2164   PetscCall(PetscFunctionListAdd(&TaoList, sname, (void (*)(void))func));
2165   PetscFunctionReturn(PETSC_SUCCESS);
2166 }
2167 
2168 /*@C
2169    TaoRegisterDestroy - Frees the list of minimization solvers that were
2170    registered by `TaoRegister()`.
2171 
2172    Not Collective
2173 
2174    Level: advanced
2175 
2176 .seealso: [](chapter_tao), `Tao`, `TaoRegisterAll()`, `TaoRegister()`
2177 @*/
2178 PetscErrorCode TaoRegisterDestroy(void)
2179 {
2180   PetscFunctionBegin;
2181   PetscCall(PetscFunctionListDestroy(&TaoList));
2182   TaoRegisterAllCalled = PETSC_FALSE;
2183   PetscFunctionReturn(PETSC_SUCCESS);
2184 }
2185 
2186 /*@
2187    TaoGetIterationNumber - Gets the number of `TaoSolve()` iterations completed
2188    at this time.
2189 
2190    Not Collective
2191 
2192    Input Parameter:
2193 .  tao - the `Tao` context
2194 
2195    Output Parameter:
2196 .  iter - iteration number
2197 
2198    Notes:
2199    For example, during the computation of iteration 2 this would return 1.
2200 
2201    Level: intermediate
2202 
2203 .seealso: [](chapter_tao), `Tao`, `TaoGetLinearSolveIterations()`, `TaoGetResidualNorm()`, `TaoGetObjective()`
2204 @*/
2205 PetscErrorCode TaoGetIterationNumber(Tao tao, PetscInt *iter)
2206 {
2207   PetscFunctionBegin;
2208   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2209   PetscValidIntPointer(iter, 2);
2210   *iter = tao->niter;
2211   PetscFunctionReturn(PETSC_SUCCESS);
2212 }
2213 
2214 /*@
2215    TaoGetResidualNorm - Gets the current value of the norm of the residual (gradient)
2216    at this time.
2217 
2218    Not Collective
2219 
2220    Input Parameter:
2221 .  tao - the `Tao` context
2222 
2223    Output Parameter:
2224 .  value - the current value
2225 
2226    Level: intermediate
2227 
2228    Developer Note:
2229    This is the 2-norm of the residual, we cannot use `TaoGetGradientNorm()` because that has
2230    a different meaning. For some reason `Tao` sometimes calls the gradient the residual.
2231 
2232 .seealso: [](chapter_tao), `Tao`, `TaoGetLinearSolveIterations()`, `TaoGetIterationNumber()`, `TaoGetObjective()`
2233 @*/
2234 PetscErrorCode TaoGetResidualNorm(Tao tao, PetscReal *value)
2235 {
2236   PetscFunctionBegin;
2237   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2238   PetscValidRealPointer(value, 2);
2239   *value = tao->residual;
2240   PetscFunctionReturn(PETSC_SUCCESS);
2241 }
2242 
2243 /*@
2244    TaoSetIterationNumber - Sets the current iteration number.
2245 
2246    Logically Collective
2247 
2248    Input Parameters:
2249 +  tao - the `Tao` context
2250 -  iter - iteration number
2251 
2252    Level: developer
2253 
2254 .seealso: [](chapter_tao), `Tao`, `TaoGetLinearSolveIterations()`
2255 @*/
2256 PetscErrorCode TaoSetIterationNumber(Tao tao, PetscInt iter)
2257 {
2258   PetscFunctionBegin;
2259   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2260   PetscValidLogicalCollectiveInt(tao, iter, 2);
2261   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2262   tao->niter = iter;
2263   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2264   PetscFunctionReturn(PETSC_SUCCESS);
2265 }
2266 
2267 /*@
2268    TaoGetTotalIterationNumber - Gets the total number of `TaoSolve()` iterations
2269    completed. This number keeps accumulating if multiple solves
2270    are called with the `Tao` object.
2271 
2272    Not Collective
2273 
2274    Input Parameter:
2275 .  tao - the `Tao` context
2276 
2277    Output Parameter:
2278 .  iter - number of iterations
2279 
2280    Level: intermediate
2281 
2282    Notes:
2283    The total iteration count is updated after each solve, if there is a current
2284    `TaoSolve()` in progress then those iterations are not included in the count
2285 
2286 .seealso: [](chapter_tao), `Tao`, `TaoGetLinearSolveIterations()`
2287 @*/
2288 PetscErrorCode TaoGetTotalIterationNumber(Tao tao, PetscInt *iter)
2289 {
2290   PetscFunctionBegin;
2291   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2292   PetscValidIntPointer(iter, 2);
2293   *iter = tao->ntotalits;
2294   PetscFunctionReturn(PETSC_SUCCESS);
2295 }
2296 
2297 /*@
2298    TaoSetTotalIterationNumber - Sets the current total iteration number.
2299 
2300    Logically Collective
2301 
2302    Input Parameters:
2303 +  tao - the `Tao` context
2304 -  iter - the iteration number
2305 
2306    Level: developer
2307 
2308 .seealso: [](chapter_tao), `Tao`, `TaoGetLinearSolveIterations()`
2309 @*/
2310 PetscErrorCode TaoSetTotalIterationNumber(Tao tao, PetscInt iter)
2311 {
2312   PetscFunctionBegin;
2313   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2314   PetscValidLogicalCollectiveInt(tao, iter, 2);
2315   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2316   tao->ntotalits = iter;
2317   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2318   PetscFunctionReturn(PETSC_SUCCESS);
2319 }
2320 
2321 /*@
2322   TaoSetConvergedReason - Sets the termination flag on a `Tao` object
2323 
2324   Logically Collective
2325 
2326   Input Parameters:
2327 + tao - the `Tao` context
2328 - reason - the `TaoConvergedReason`
2329 
2330    Level: intermediate
2331 
2332 .seealso: [](chapter_tao), `Tao`, `TaoConvergedReason`
2333 @*/
2334 PetscErrorCode TaoSetConvergedReason(Tao tao, TaoConvergedReason reason)
2335 {
2336   PetscFunctionBegin;
2337   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2338   PetscValidLogicalCollectiveEnum(tao, reason, 2);
2339   tao->reason = reason;
2340   PetscFunctionReturn(PETSC_SUCCESS);
2341 }
2342 
2343 /*@
2344    TaoGetConvergedReason - Gets the reason the `TaoSolve()` was stopped.
2345 
2346    Not Collective
2347 
2348    Input Parameter:
2349 .  tao - the `Tao` solver context
2350 
2351    Output Parameter:
2352 .  reason - value of `TaoConvergedReason`
2353 
2354    Level: intermediate
2355 
2356 .seealso: [](chapter_tao), `Tao`, `TaoConvergedReason`, `TaoSetConvergenceTest()`, `TaoSetTolerances()`
2357 @*/
2358 PetscErrorCode TaoGetConvergedReason(Tao tao, TaoConvergedReason *reason)
2359 {
2360   PetscFunctionBegin;
2361   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2362   PetscValidPointer(reason, 2);
2363   *reason = tao->reason;
2364   PetscFunctionReturn(PETSC_SUCCESS);
2365 }
2366 
2367 /*@
2368    TaoGetSolutionStatus - Get the current iterate, objective value,
2369    residual, infeasibility, and termination from a `Tao` object
2370 
2371    Not Collective
2372 
2373    Input Parameter:
2374 .  tao - the `Tao` context
2375 
2376    Output Parameters:
2377 +  iterate - the current iterate number (>=0)
2378 .  f - the current function value
2379 .  gnorm - the square of the gradient norm, duality gap, or other measure indicating distance from optimality.
2380 .  cnorm - the infeasibility of the current solution with regard to the constraints.
2381 .  xdiff - the step length or trust region radius of the most recent iterate.
2382 -  reason - The termination reason, which can equal `TAO_CONTINUE_ITERATING`
2383 
2384    Level: intermediate
2385 
2386    Notes:
2387    Tao returns the values set by the solvers in the routine `TaoMonitor()`.
2388 
2389    If any of the output arguments are set to `NULL`, no corresponding value will be returned.
2390 
2391 .seealso: [](chapter_tao), `TaoMonitor()`, `TaoGetConvergedReason()`
2392 @*/
2393 PetscErrorCode TaoGetSolutionStatus(Tao tao, PetscInt *its, PetscReal *f, PetscReal *gnorm, PetscReal *cnorm, PetscReal *xdiff, TaoConvergedReason *reason)
2394 {
2395   PetscFunctionBegin;
2396   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2397   if (its) *its = tao->niter;
2398   if (f) *f = tao->fc;
2399   if (gnorm) *gnorm = tao->residual;
2400   if (cnorm) *cnorm = tao->cnorm;
2401   if (reason) *reason = tao->reason;
2402   if (xdiff) *xdiff = tao->step;
2403   PetscFunctionReturn(PETSC_SUCCESS);
2404 }
2405 
2406 /*@C
2407    TaoGetType - Gets the current `TaoType` being used in the `Tao` object
2408 
2409    Not Collective
2410 
2411    Input Parameter:
2412 .  tao - the `Tao` solver context
2413 
2414    Output Parameter:
2415 .  type - the `TaoType`
2416 
2417    Level: intermediate
2418 
2419 .seealso: [](chapter_tao), `Tao`, `TaoType`, `TaoSetType()`
2420 @*/
2421 PetscErrorCode TaoGetType(Tao tao, TaoType *type)
2422 {
2423   PetscFunctionBegin;
2424   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2425   PetscValidPointer(type, 2);
2426   *type = ((PetscObject)tao)->type_name;
2427   PetscFunctionReturn(PETSC_SUCCESS);
2428 }
2429 
2430 /*@C
2431   TaoMonitor - Monitor the solver and the current solution.  This
2432   routine will record the iteration number and residual statistics,
2433   and call any monitors specified by the user.
2434 
2435    Input Parameters:
2436 +  tao - the `Tao` context
2437 .  its - the current iterate number (>=0)
2438 .  f - the current objective function value
2439 .  res - the gradient norm, square root of the duality gap, or other measure indicating distince from optimality.  This measure will be recorded and
2440           used for some termination tests.
2441 .  cnorm - the infeasibility of the current solution with regard to the constraints.
2442 -  steplength - multiple of the step direction added to the previous iterate.
2443 
2444    Output Parameter:
2445 .  reason - The termination reason, which can equal `TAO_CONTINUE_ITERATING`
2446 
2447    Options Database Key:
2448 .  -tao_monitor - Use the default monitor, which prints statistics to standard output
2449 
2450    Level: developer
2451 
2452 .seealso: [](chapter_tao), `Tao`, `TaoGetConvergedReason()`, `TaoMonitorDefault()`, `TaoSetMonitor()`
2453 @*/
2454 PetscErrorCode TaoMonitor(Tao tao, PetscInt its, PetscReal f, PetscReal res, PetscReal cnorm, PetscReal steplength)
2455 {
2456   PetscInt i;
2457 
2458   PetscFunctionBegin;
2459   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2460   tao->fc       = f;
2461   tao->residual = res;
2462   tao->cnorm    = cnorm;
2463   tao->step     = steplength;
2464   if (!its) {
2465     tao->cnorm0 = cnorm;
2466     tao->gnorm0 = res;
2467   }
2468   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(res), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
2469   for (i = 0; i < tao->numbermonitors; i++) PetscCall((*tao->monitor[i])(tao, tao->monitorcontext[i]));
2470   PetscFunctionReturn(PETSC_SUCCESS);
2471 }
2472 
2473 /*@
2474    TaoSetConvergenceHistory - Sets the array used to hold the convergence history.
2475 
2476    Logically Collective
2477 
2478    Input Parameters:
2479 +  tao - the `Tao` solver context
2480 .  obj   - array to hold objective value history
2481 .  resid - array to hold residual history
2482 .  cnorm - array to hold constraint violation history
2483 .  lits - integer array holds the number of linear iterations for each Tao iteration
2484 .  na  - size of `obj`, `resid`, and `cnorm`
2485 -  reset - `PETSC_TRUE` indicates each new minimization resets the history counter to zero,
2486            else it continues storing new values for new minimizations after the old ones
2487 
2488    Level: intermediate
2489 
2490    Notes:
2491    If set, `Tao` will fill the given arrays with the indicated
2492    information at each iteration.  If 'obj','resid','cnorm','lits' are
2493    *all* `NULL` then space (using size `na`, or 1000 if na is `PETSC_DECIDE` or
2494    `PETSC_DEFAULT`) is allocated for the history.
2495    If not all are `NULL`, then only the non-`NULL` information categories
2496    will be stored, the others will be ignored.
2497 
2498    Any convergence information after iteration number 'na' will not be stored.
2499 
2500    This routine is useful, e.g., when running a code for purposes
2501    of accurate performance monitoring, when no I/O should be done
2502    during the section of code that is being timed.
2503 
2504 .seealso: [](chapter_tao), `TaoGetConvergenceHistory()`
2505 @*/
2506 PetscErrorCode TaoSetConvergenceHistory(Tao tao, PetscReal obj[], PetscReal resid[], PetscReal cnorm[], PetscInt lits[], PetscInt na, PetscBool reset)
2507 {
2508   PetscFunctionBegin;
2509   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2510   if (obj) PetscValidRealPointer(obj, 2);
2511   if (resid) PetscValidRealPointer(resid, 3);
2512   if (cnorm) PetscValidRealPointer(cnorm, 4);
2513   if (lits) PetscValidIntPointer(lits, 5);
2514 
2515   if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
2516   if (!obj && !resid && !cnorm && !lits) {
2517     PetscCall(PetscCalloc4(na, &obj, na, &resid, na, &cnorm, na, &lits));
2518     tao->hist_malloc = PETSC_TRUE;
2519   }
2520 
2521   tao->hist_obj   = obj;
2522   tao->hist_resid = resid;
2523   tao->hist_cnorm = cnorm;
2524   tao->hist_lits  = lits;
2525   tao->hist_max   = na;
2526   tao->hist_reset = reset;
2527   tao->hist_len   = 0;
2528   PetscFunctionReturn(PETSC_SUCCESS);
2529 }
2530 
2531 /*@C
2532    TaoGetConvergenceHistory - Gets the arrays used that hold the convergence history.
2533 
2534    Collective
2535 
2536    Input Parameter:
2537 .  tao - the `Tao` context
2538 
2539    Output Parameters:
2540 +  obj   - array used to hold objective value history
2541 .  resid - array used to hold residual history
2542 .  cnorm - array used to hold constraint violation history
2543 .  lits  - integer array used to hold linear solver iteration count
2544 -  nhist  - size of `obj`, `resid`, `cnorm`, and `lits`
2545 
2546    Level: advanced
2547 
2548    Notes:
2549     This routine must be preceded by calls to `TaoSetConvergenceHistory()`
2550     and `TaoSolve()`, otherwise it returns useless information.
2551 
2552    This routine is useful, e.g., when running a code for purposes
2553    of accurate performance monitoring, when no I/O should be done
2554    during the section of code that is being timed.
2555 
2556    Fortran Note:
2557    The calling sequence is
2558 .vb
2559    call TaoGetConvergenceHistory(Tao tao, PetscInt nhist, PetscErrorCode ierr)
2560 .ve
2561 
2562 .seealso: [](chapter_tao), `Tao`, `TaoSolve()`, `TaoSetConvergenceHistory()`
2563 @*/
2564 PetscErrorCode TaoGetConvergenceHistory(Tao tao, PetscReal **obj, PetscReal **resid, PetscReal **cnorm, PetscInt **lits, PetscInt *nhist)
2565 {
2566   PetscFunctionBegin;
2567   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2568   if (obj) *obj = tao->hist_obj;
2569   if (cnorm) *cnorm = tao->hist_cnorm;
2570   if (resid) *resid = tao->hist_resid;
2571   if (lits) *lits = tao->hist_lits;
2572   if (nhist) *nhist = tao->hist_len;
2573   PetscFunctionReturn(PETSC_SUCCESS);
2574 }
2575 
2576 /*@
2577    TaoSetApplicationContext - Sets the optional user-defined context for a `Tao` solver.
2578 
2579    Logically Collective
2580 
2581    Input Parameters:
2582 +  tao  - the `Tao` context
2583 -  usrP - optional user context
2584 
2585    Level: intermediate
2586 
2587 .seealso: [](chapter_tao), `Tao`, `TaoGetApplicationContext()`, `TaoSetApplicationContext()`
2588 @*/
2589 PetscErrorCode TaoSetApplicationContext(Tao tao, void *usrP)
2590 {
2591   PetscFunctionBegin;
2592   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2593   tao->user = usrP;
2594   PetscFunctionReturn(PETSC_SUCCESS);
2595 }
2596 
2597 /*@
2598    TaoGetApplicationContext - Gets the user-defined context for a `Tao` solver
2599 
2600    Not Collective
2601 
2602    Input Parameter:
2603 .  tao  - the `Tao` context
2604 
2605    Output Parameter:
2606 .  usrP - user context
2607 
2608    Level: intermediate
2609 
2610 .seealso: [](chapter_tao), `Tao`, `TaoSetApplicationContext()`
2611 @*/
2612 PetscErrorCode TaoGetApplicationContext(Tao tao, void *usrP)
2613 {
2614   PetscFunctionBegin;
2615   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2616   PetscValidPointer(usrP, 2);
2617   *(void **)usrP = tao->user;
2618   PetscFunctionReturn(PETSC_SUCCESS);
2619 }
2620 
2621 /*@
2622    TaoSetGradientNorm - Sets the matrix used to define the norm that measures the size of the gradient.
2623 
2624    Collective
2625 
2626    Input Parameters:
2627 +  tao  - the `Tao` context
2628 -  M    - matrix that defines the norm
2629 
2630    Level: beginner
2631 
2632 .seealso: [](chapter_tao), `Tao`, `TaoGetGradientNorm()`, `TaoGradientNorm()`
2633 @*/
2634 PetscErrorCode TaoSetGradientNorm(Tao tao, Mat M)
2635 {
2636   PetscFunctionBegin;
2637   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2638   PetscValidHeaderSpecific(M, MAT_CLASSID, 2);
2639   PetscCall(PetscObjectReference((PetscObject)M));
2640   PetscCall(MatDestroy(&tao->gradient_norm));
2641   PetscCall(VecDestroy(&tao->gradient_norm_tmp));
2642   tao->gradient_norm = M;
2643   PetscCall(MatCreateVecs(M, NULL, &tao->gradient_norm_tmp));
2644   PetscFunctionReturn(PETSC_SUCCESS);
2645 }
2646 
2647 /*@
2648    TaoGetGradientNorm - Returns the matrix used to define the norm used for measuring the size of the gradient.
2649 
2650    Not Collective
2651 
2652    Input Parameter:
2653 .  tao  - the `Tao` context
2654 
2655    Output Parameter:
2656 .  M - gradient norm
2657 
2658    Level: beginner
2659 
2660 .seealso: [](chapter_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGradientNorm()`
2661 @*/
2662 PetscErrorCode TaoGetGradientNorm(Tao tao, Mat *M)
2663 {
2664   PetscFunctionBegin;
2665   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2666   PetscValidPointer(M, 2);
2667   *M = tao->gradient_norm;
2668   PetscFunctionReturn(PETSC_SUCCESS);
2669 }
2670 
2671 /*@C
2672    TaoGradientNorm - Compute the norm using the `NormType`, the user has selected
2673 
2674    Collective
2675 
2676    Input Parameters:
2677 +  tao      - the `Tao` context
2678 .  gradient - the gradient to be computed
2679 -  norm     - the norm type
2680 
2681    Output Parameter:
2682 .  gnorm    - the gradient norm
2683 
2684    Level: advanced
2685 
2686 .seealso: [](chapter_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGetGradientNorm()`
2687 @*/
2688 PetscErrorCode TaoGradientNorm(Tao tao, Vec gradient, NormType type, PetscReal *gnorm)
2689 {
2690   PetscFunctionBegin;
2691   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2692   PetscValidHeaderSpecific(gradient, VEC_CLASSID, 2);
2693   PetscValidLogicalCollectiveEnum(tao, type, 3);
2694   PetscValidRealPointer(gnorm, 4);
2695   if (tao->gradient_norm) {
2696     PetscScalar gnorms;
2697 
2698     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.");
2699     PetscCall(MatMult(tao->gradient_norm, gradient, tao->gradient_norm_tmp));
2700     PetscCall(VecDot(gradient, tao->gradient_norm_tmp, &gnorms));
2701     *gnorm = PetscRealPart(PetscSqrtScalar(gnorms));
2702   } else {
2703     PetscCall(VecNorm(gradient, type, gnorm));
2704   }
2705   PetscFunctionReturn(PETSC_SUCCESS);
2706 }
2707 
2708 /*@C
2709    TaoMonitorDrawCtxCreate - Creates the monitor context for `TaoMonitorDrawSolution()`
2710 
2711    Collective
2712 
2713    Output Parameter:
2714 .    ctx - the monitor context
2715 
2716    Options Database Key:
2717 .   -tao_draw_solution_initial - show initial guess as well as current solution
2718 
2719    Level: intermediate
2720 
2721 .seealso: [](chapter_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorDrawCtx()`
2722 @*/
2723 PetscErrorCode TaoMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TaoMonitorDrawCtx *ctx)
2724 {
2725   PetscFunctionBegin;
2726   PetscCall(PetscNew(ctx));
2727   PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
2728   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
2729   (*ctx)->howoften = howoften;
2730   PetscFunctionReturn(PETSC_SUCCESS);
2731 }
2732 
2733 /*@C
2734    TaoMonitorDrawCtxDestroy - Destroys the monitor context for `TaoMonitorDrawSolution()`
2735 
2736    Collective
2737 
2738    Input Parameter:
2739 .    ctx - the monitor context
2740 
2741    Level: intermediate
2742 
2743 .seealso: [](chapter_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorDrawSolution()`
2744 @*/
2745 PetscErrorCode TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx *ictx)
2746 {
2747   PetscFunctionBegin;
2748   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
2749   PetscCall(PetscFree(*ictx));
2750   PetscFunctionReturn(PETSC_SUCCESS);
2751 }
2752