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