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