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