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