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