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