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