Actual source code: taosolver.c

  1: #include <petsc/private/taoimpl.h>
  2: #include <petsc/private/snesimpl.h>

  4: PetscBool         TaoRegisterAllCalled = PETSC_FALSE;
  5: PetscFunctionList TaoList              = NULL;

  7: PetscClassId TAO_CLASSID;

  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;

 17: const char *TaoSubSetTypes[] = {"subvec", "mask", "matrixfree", "TaoSubSetType", "TAO_SUBSET_", NULL};

 19: struct _n_TaoMonitorDrawCtx {
 20:   PetscViewer viewer;
 21:   PetscInt    howoften; /* when > 0 uses iteration % howoften, when negative only final solution plotted */
 22: };

 24: static PetscErrorCode KSPPreSolve_TAOEW_Private(KSP ksp, Vec b, Vec x, PetscCtx ctx)
 25: {
 26:   Tao  tao          = (Tao)ctx;
 27:   SNES snes_ewdummy = tao->snes_ewdummy;

 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: }

 41: static PetscErrorCode KSPPostSolve_TAOEW_Private(KSP ksp, Vec b, Vec x, PetscCtx ctx)
 42: {
 43:   Tao  tao          = (Tao)ctx;
 44:   SNES snes_ewdummy = tao->snes_ewdummy;

 46:   PetscFunctionBegin;
 47:   if (!snes_ewdummy) PetscFunctionReturn(PETSC_SUCCESS);
 48:   PetscCall(KSPPostSolve_SNESEW(ksp, b, x, snes_ewdummy));
 49:   PetscFunctionReturn(PETSC_SUCCESS);
 50: }

 52: static PetscErrorCode TaoSetUpEW_Private(Tao tao)
 53: {
 54:   SNESKSPEW  *kctx;
 55:   const char *ewprefix;

 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));

 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: }

 72: /*@
 73:   TaoParametersInitialize - Sets all the parameters in `tao` to their default value (when `TaoCreate()` was called) if they
 74:   currently contain default values. Default values are the parameter values when the object's type is set.

 76:   Collective

 78:   Input Parameter:
 79: . tao - the `Tao` object

 81:   Level: developer

 83:   Developer Note:
 84:   This is called by all the `TaoCreate_XXX()` routines.

 86: .seealso: [](ch_snes), `Tao`, `TaoSolve()`, `TaoDestroy()`,
 87:           `PetscObjectParameterSetDefault()`
 88: @*/
 89: PetscErrorCode TaoParametersInitialize(Tao tao)
 90: {
 91:   PetscObjectParameterSetDefault(tao, max_it, 10000);
 92:   PetscObjectParameterSetDefault(tao, max_funcs, PETSC_UNLIMITED);
 93:   PetscObjectParameterSetDefault(tao, gatol, PetscDefined(USE_REAL_SINGLE) ? 1e-5 : 1e-8);
 94:   PetscObjectParameterSetDefault(tao, grtol, PetscDefined(USE_REAL_SINGLE) ? 1e-5 : 1e-8);
 95:   PetscObjectParameterSetDefault(tao, crtol, PetscDefined(USE_REAL_SINGLE) ? 1e-5 : 1e-8);
 96:   PetscObjectParameterSetDefault(tao, catol, PetscDefined(USE_REAL_SINGLE) ? 1e-5 : 1e-8);
 97:   PetscObjectParameterSetDefault(tao, gttol, 0.0);
 98:   PetscObjectParameterSetDefault(tao, steptol, 0.0);
 99:   PetscObjectParameterSetDefault(tao, fmin, PETSC_NINFINITY);
100:   PetscObjectParameterSetDefault(tao, trust0, PETSC_INFINITY);
101:   return PETSC_SUCCESS;
102: }

104: /*@
105:   TaoCreate - Creates a Tao solver

107:   Collective

109:   Input Parameter:
110: . comm - MPI communicator

112:   Output Parameter:
113: . newtao - the new `Tao` context

115:   Options Database Key:
116: . -tao_type - select which method Tao should use

118:   Level: beginner

120: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoDestroy()`, `TaoSetFromOptions()`, `TaoSetType()`
121: @*/
122: PetscErrorCode TaoCreate(MPI_Comm comm, Tao *newtao)
123: {
124:   Tao tao;

126:   PetscFunctionBegin;
127:   PetscAssertPointer(newtao, 2);
128:   PetscCall(TaoInitializePackage());
129:   PetscCall(TaoLineSearchInitializePackage());

131:   PetscCall(PetscHeaderCreate(tao, TAO_CLASSID, "Tao", "Optimization solver", "Tao", comm, TaoDestroy, TaoView));
132:   tao->ops->convergencetest = TaoDefaultConvergenceTest;

134:   tao->hist_reset = PETSC_TRUE;
135:   PetscCall(TaoResetStatistics(tao));
136:   *newtao = tao;
137:   PetscFunctionReturn(PETSC_SUCCESS);
138: }

140: /*@
141:   TaoSolve - Solves an optimization problem min F(x) s.t. l <= x <= u

143:   Collective

145:   Input Parameter:
146: . tao - the `Tao` context

148:   Level: beginner

150:   Notes:
151:   The user must set up the `Tao` object  with calls to `TaoSetSolution()`, `TaoSetObjective()`, `TaoSetGradient()`, and (if using 2nd order method) `TaoSetHessian()`.

153:   You should call `TaoGetConvergedReason()` or run with `-tao_converged_reason` to determine if the optimization algorithm actually succeeded or
154:   why it failed.

156: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSetObjective()`, `TaoSetGradient()`, `TaoSetHessian()`, `TaoGetConvergedReason()`, `TaoSetUp()`
157:  @*/
158: PetscErrorCode TaoSolve(Tao tao)
159: {
160:   static PetscBool set = PETSC_FALSE;

162:   PetscFunctionBegin;
164:   PetscCall(PetscCitationsRegister("@TechReport{tao-user-ref,\n"
165:                                    "title   = {Toolkit for Advanced Optimization (TAO) Users Manual},\n"
166:                                    "author  = {Todd Munson and Jason Sarich and Stefan Wild and Steve Benson and Lois Curfman McInnes},\n"
167:                                    "Institution = {Argonne National Laboratory},\n"
168:                                    "Year   = 2014,\n"
169:                                    "Number = {ANL/MCS-TM-322 - Revision 3.5},\n"
170:                                    "url    = {https://www.mcs.anl.gov/research/projects/tao/}\n}\n",
171:                                    &set));
172:   tao->header_printed = PETSC_FALSE;
173:   PetscCall(TaoSetUp(tao));
174:   PetscCall(TaoResetStatistics(tao));
175:   if (tao->linesearch) PetscCall(TaoLineSearchReset(tao->linesearch));

177:   PetscCall(PetscLogEventBegin(TAO_Solve, tao, 0, 0, 0));
178:   PetscTryTypeMethod(tao, solve);
179:   PetscCall(PetscLogEventEnd(TAO_Solve, tao, 0, 0, 0));

181:   PetscCall(VecViewFromOptions(tao->solution, (PetscObject)tao, "-tao_view_solution"));

183:   tao->ntotalits += tao->niter;

185:   if (tao->printreason) {
186:     PetscViewer viewer = PETSC_VIEWER_STDOUT_(((PetscObject)tao)->comm);

188:     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)tao)->tablevel));
189:     if (tao->reason > 0) {
190:       if (((PetscObject)tao)->prefix) {
191:         PetscCall(PetscViewerASCIIPrintf(viewer, "  TAO %s solve converged due to %s iterations %" PetscInt_FMT "\n", ((PetscObject)tao)->prefix, TaoConvergedReasons[tao->reason], tao->niter));
192:       } else {
193:         PetscCall(PetscViewerASCIIPrintf(viewer, "  TAO solve converged due to %s iterations %" PetscInt_FMT "\n", TaoConvergedReasons[tao->reason], tao->niter));
194:       }
195:     } else {
196:       if (((PetscObject)tao)->prefix) {
197:         PetscCall(PetscViewerASCIIPrintf(viewer, "  TAO %s solve did not converge due to %s iteration %" PetscInt_FMT "\n", ((PetscObject)tao)->prefix, TaoConvergedReasons[tao->reason], tao->niter));
198:       } else {
199:         PetscCall(PetscViewerASCIIPrintf(viewer, "  TAO solve did not converge due to %s iteration %" PetscInt_FMT "\n", TaoConvergedReasons[tao->reason], tao->niter));
200:       }
201:     }
202:     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)tao)->tablevel));
203:   }
204:   PetscCall(TaoViewFromOptions(tao, NULL, "-tao_view"));
205:   PetscFunctionReturn(PETSC_SUCCESS);
206: }

208: /*@
209:   TaoSetUp - Sets up the internal data structures for the later use
210:   of a Tao solver

212:   Collective

214:   Input Parameter:
215: . tao - the `Tao` context

217:   Level: advanced

219:   Note:
220:   The user will not need to explicitly call `TaoSetUp()`, as it will
221:   automatically be called in `TaoSolve()`.  However, if the user
222:   desires to call it explicitly, it should come after `TaoCreate()`
223:   and any TaoSetSomething() routines, but before `TaoSolve()`.

225: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
226: @*/
227: PetscErrorCode TaoSetUp(Tao tao)
228: {
229:   PetscFunctionBegin;
231:   if (tao->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
232:   PetscCall(TaoSetUpEW_Private(tao));
233:   PetscCheck(tao->solution, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoSetSolution");
234:   PetscTryTypeMethod(tao, setup);
235:   tao->setupcalled = PETSC_TRUE;
236:   PetscFunctionReturn(PETSC_SUCCESS);
237: }

239: /*@
240:   TaoDestroy - Destroys the `Tao` context that was created with `TaoCreate()`

242:   Collective

244:   Input Parameter:
245: . tao - the `Tao` context

247:   Level: beginner

249: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
250: @*/
251: PetscErrorCode TaoDestroy(Tao *tao)
252: {
253:   PetscFunctionBegin;
254:   if (!*tao) PetscFunctionReturn(PETSC_SUCCESS);
256:   if (--((PetscObject)*tao)->refct > 0) {
257:     *tao = NULL;
258:     PetscFunctionReturn(PETSC_SUCCESS);
259:   }

261:   PetscTryTypeMethod(*tao, destroy);
262:   PetscCall(KSPDestroy(&(*tao)->ksp));
263:   PetscCall(SNESDestroy(&(*tao)->snes_ewdummy));
264:   PetscCall(TaoLineSearchDestroy(&(*tao)->linesearch));

266:   if ((*tao)->ops->convergencedestroy) {
267:     PetscCall((*(*tao)->ops->convergencedestroy)((*tao)->cnvP));
268:     if ((*tao)->jacobian_state_inv) PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
269:   }
270:   PetscCall(VecDestroy(&(*tao)->solution));
271:   PetscCall(VecDestroy(&(*tao)->gradient));
272:   PetscCall(VecDestroy(&(*tao)->ls_res));

274:   if ((*tao)->gradient_norm) {
275:     PetscCall(PetscObjectDereference((PetscObject)(*tao)->gradient_norm));
276:     PetscCall(VecDestroy(&(*tao)->gradient_norm_tmp));
277:   }

279:   PetscCall(VecDestroy(&(*tao)->XL));
280:   PetscCall(VecDestroy(&(*tao)->XU));
281:   PetscCall(VecDestroy(&(*tao)->IL));
282:   PetscCall(VecDestroy(&(*tao)->IU));
283:   PetscCall(VecDestroy(&(*tao)->DE));
284:   PetscCall(VecDestroy(&(*tao)->DI));
285:   PetscCall(VecDestroy(&(*tao)->constraints));
286:   PetscCall(VecDestroy(&(*tao)->constraints_equality));
287:   PetscCall(VecDestroy(&(*tao)->constraints_inequality));
288:   PetscCall(VecDestroy(&(*tao)->stepdirection));
289:   PetscCall(MatDestroy(&(*tao)->hessian_pre));
290:   PetscCall(MatDestroy(&(*tao)->hessian));
291:   PetscCall(MatDestroy(&(*tao)->ls_jac));
292:   PetscCall(MatDestroy(&(*tao)->ls_jac_pre));
293:   PetscCall(MatDestroy(&(*tao)->jacobian_pre));
294:   PetscCall(MatDestroy(&(*tao)->jacobian));
295:   PetscCall(MatDestroy(&(*tao)->jacobian_state_pre));
296:   PetscCall(MatDestroy(&(*tao)->jacobian_state));
297:   PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
298:   PetscCall(MatDestroy(&(*tao)->jacobian_design));
299:   PetscCall(MatDestroy(&(*tao)->jacobian_equality));
300:   PetscCall(MatDestroy(&(*tao)->jacobian_equality_pre));
301:   PetscCall(MatDestroy(&(*tao)->jacobian_inequality));
302:   PetscCall(MatDestroy(&(*tao)->jacobian_inequality_pre));
303:   PetscCall(ISDestroy(&(*tao)->state_is));
304:   PetscCall(ISDestroy(&(*tao)->design_is));
305:   PetscCall(VecDestroy(&(*tao)->res_weights_v));
306:   PetscCall(TaoMonitorCancel(*tao));
307:   if ((*tao)->hist_malloc) PetscCall(PetscFree4((*tao)->hist_obj, (*tao)->hist_resid, (*tao)->hist_cnorm, (*tao)->hist_lits));
308:   if ((*tao)->res_weights_n) {
309:     PetscCall(PetscFree((*tao)->res_weights_rows));
310:     PetscCall(PetscFree((*tao)->res_weights_cols));
311:     PetscCall(PetscFree((*tao)->res_weights_w));
312:   }
313:   PetscCall(PetscHeaderDestroy(tao));
314:   PetscFunctionReturn(PETSC_SUCCESS);
315: }

317: /*@
318:   TaoKSPSetUseEW - Sets `SNES` to use Eisenstat-Walker method {cite}`ew96`for computing relative tolerance for linear solvers.

320:   Logically Collective

322:   Input Parameters:
323: + tao  - Tao context
324: - flag - `PETSC_TRUE` or `PETSC_FALSE`

326:   Level: advanced

328:   Note:
329:   See `SNESKSPSetUseEW()` for customization details.

331: .seealso: [](ch_tao), `Tao`, `SNESKSPSetUseEW()`
332: @*/
333: PetscErrorCode TaoKSPSetUseEW(Tao tao, PetscBool flag)
334: {
335:   PetscFunctionBegin;
338:   tao->ksp_ewconv = flag;
339:   PetscFunctionReturn(PETSC_SUCCESS);
340: }

342: /*@C
343:   TaoMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user

345:   Collective

347:   Input Parameters:
348: + tao     - `Tao` object you wish to monitor
349: . name    - the monitor type one is seeking
350: . help    - message indicating what monitoring is done
351: . manual  - manual page for the monitor
352: - monitor - the monitor function, this must use a `PetscViewerFormat` as its context

354:   Level: developer

356: .seealso: [](ch_ts), `Tao`, `TaoMonitorSet()`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
357:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`
358:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
359:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
360:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
361:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
362:           `PetscOptionsFList()`, `PetscOptionsEList()`
363: @*/
364: PetscErrorCode TaoMonitorSetFromOptions(Tao tao, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(Tao, PetscViewerAndFormat *))
365: {
366:   PetscViewer       viewer;
367:   PetscViewerFormat format;
368:   PetscBool         flg;

370:   PetscFunctionBegin;
371:   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)tao), ((PetscObject)tao)->options, ((PetscObject)tao)->prefix, name, &viewer, &format, &flg));
372:   if (flg) {
373:     PetscViewerAndFormat *vf;
374:     char                  interval_key[1024];

376:     PetscCall(PetscSNPrintf(interval_key, sizeof interval_key, "%s_interval", name));
377:     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
378:     vf->view_interval = 1;
379:     PetscCall(PetscOptionsGetInt(((PetscObject)tao)->options, ((PetscObject)tao)->prefix, interval_key, &vf->view_interval, NULL));

381:     PetscCall(PetscViewerDestroy(&viewer));
382:     PetscCall(TaoMonitorSet(tao, (PetscErrorCode (*)(Tao, PetscCtx))monitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy));
383:   }
384:   PetscFunctionReturn(PETSC_SUCCESS);
385: }

387: /*@
388:   TaoSetFromOptions - Sets various Tao parameters from the options database

390:   Collective

392:   Input Parameter:
393: . tao - the `Tao` solver context

395:   Options Database Keys:
396: + -tao_type <type>             - The algorithm that Tao uses (lmvm, nls, etc.)
397: . -tao_gatol <gatol>           - absolute error tolerance for ||gradient||
398: . -tao_grtol <grtol>           - relative error tolerance for ||gradient||
399: . -tao_gttol <gttol>           - reduction of ||gradient|| relative to initial gradient
400: . -tao_max_it <max>            - sets maximum number of iterations
401: . -tao_max_funcs <max>         - sets maximum number of function evaluations
402: . -tao_fmin <fmin>             - stop if function value reaches fmin
403: . -tao_steptol <tol>           - stop if trust region radius less than <tol>
404: . -tao_trust0 <t>              - initial trust region radius
405: . -tao_view_solution           - view the solution at the end of the optimization process
406: . -tao_monitor                 - prints function value and residual norm at each iteration
407: . -tao_monitor_short           - same as `-tao_monitor`, but truncates very small values
408: . -tao_monitor_constraint_norm - prints objective value, gradient, and constraint norm at each iteration
409: . -tao_monitor_globalization   - prints information about the globalization at each iteration
410: . -tao_monitor_solution        - prints solution vector at each iteration
411: . -tao_monitor_ls_residual     - prints least-squares residual vector at each iteration
412: . -tao_monitor_step            - prints step vector at each iteration
413: . -tao_monitor_gradient        - prints gradient vector at each iteration
414: . -tao_monitor_solution_draw   - graphically view solution vector at each iteration
415: . -tao_monitor_step_draw       - graphically view step vector at each iteration
416: . -tao_monitor_gradient_draw   - graphically view gradient at each iteration
417: . -tao_monitor_cancel          - cancels all monitors (except those set with command line)
418: . -tao_fd_gradient             - use gradient computed with finite differences
419: . -tao_fd_hessian              - use hessian computed with finite differences
420: . -tao_mf_hessian              - use matrix-free Hessian computed with finite differences
421: . -tao_view                    - prints information about the Tao after solving
422: - -tao_converged_reason        - prints the reason Tao stopped iterating

424:   Level: beginner

426:   Note:
427:   To see all options, run your program with the `-help` option or consult the
428:   user's manual. Should be called after `TaoCreate()` but before `TaoSolve()`

430: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
431: @*/
432: PetscErrorCode TaoSetFromOptions(Tao tao)
433: {
434:   TaoType   default_type = TAOLMVM;
435:   char      type[256];
436:   PetscBool flg, found;
437:   MPI_Comm  comm;
438:   PetscReal catol, crtol, gatol, grtol, gttol;

440:   PetscFunctionBegin;
442:   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));

444:   if (((PetscObject)tao)->type_name) default_type = ((PetscObject)tao)->type_name;

446:   PetscObjectOptionsBegin((PetscObject)tao);
447:   /* Check for type from options */
448:   PetscCall(PetscOptionsFList("-tao_type", "Tao Solver type", "TaoSetType", TaoList, default_type, type, 256, &flg));
449:   if (flg) {
450:     PetscCall(TaoSetType(tao, type));
451:   } else if (!((PetscObject)tao)->type_name) {
452:     PetscCall(TaoSetType(tao, default_type));
453:   }

455:   /* Tao solvers do not set the prefix, set it here if not yet done
456:      We do it after SetType since solver may have been changed */
457:   if (tao->linesearch) {
458:     const char *prefix;
459:     PetscCall(TaoLineSearchGetOptionsPrefix(tao->linesearch, &prefix));
460:     if (!prefix) PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, ((PetscObject)tao)->prefix));
461:   }

463:   catol = tao->catol;
464:   crtol = tao->crtol;
465:   PetscCall(PetscOptionsReal("-tao_catol", "Stop if constraints violations within", "TaoSetConstraintTolerances", tao->catol, &catol, NULL));
466:   PetscCall(PetscOptionsReal("-tao_crtol", "Stop if relative constraint violations within", "TaoSetConstraintTolerances", tao->crtol, &crtol, NULL));
467:   PetscCall(TaoSetConstraintTolerances(tao, catol, crtol));

469:   gatol = tao->gatol;
470:   grtol = tao->grtol;
471:   gttol = tao->gttol;
472:   PetscCall(PetscOptionsReal("-tao_gatol", "Stop if norm of gradient less than", "TaoSetTolerances", tao->gatol, &gatol, NULL));
473:   PetscCall(PetscOptionsReal("-tao_grtol", "Stop if norm of gradient divided by the function value is less than", "TaoSetTolerances", tao->grtol, &grtol, NULL));
474:   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, &gttol, NULL));
475:   PetscCall(TaoSetTolerances(tao, gatol, grtol, gttol));

477:   PetscCall(PetscOptionsInt("-tao_max_it", "Stop if iteration number exceeds", "TaoSetMaximumIterations", tao->max_it, &tao->max_it, &flg));
478:   if (flg) PetscCall(TaoSetMaximumIterations(tao, tao->max_it));

480:   PetscCall(PetscOptionsInt("-tao_max_funcs", "Stop if number of function evaluations exceeds", "TaoSetMaximumFunctionEvaluations", tao->max_funcs, &tao->max_funcs, &flg));
481:   if (flg) PetscCall(TaoSetMaximumFunctionEvaluations(tao, tao->max_funcs));

483:   PetscCall(PetscOptionsReal("-tao_fmin", "Stop if function less than", "TaoSetFunctionLowerBound", tao->fmin, &tao->fmin, NULL));
484:   PetscCall(PetscOptionsBoundedReal("-tao_steptol", "Stop if step size or trust region radius less than", "", tao->steptol, &tao->steptol, NULL, 0));
485:   PetscCall(PetscOptionsReal("-tao_trust0", "Initial trust region radius", "TaoSetInitialTrustRegionRadius", tao->trust0, &tao->trust0, &flg));
486:   if (flg) PetscCall(TaoSetInitialTrustRegionRadius(tao, tao->trust0));

488:   PetscCall(PetscOptionsDeprecated("-tao_solution_monitor", "-tao_monitor_solution", "3.21", NULL));
489:   PetscCall(PetscOptionsDeprecated("-tao_gradient_monitor", "-tao_monitor_gradient", "3.21", NULL));
490:   PetscCall(PetscOptionsDeprecated("-tao_stepdirection_monitor", "-tao_monitor_step", "3.21", NULL));
491:   PetscCall(PetscOptionsDeprecated("-tao_residual_monitor", "-tao_monitor_residual", "3.21", NULL));
492:   PetscCall(PetscOptionsDeprecated("-tao_smonitor", "-tao_monitor_short", "3.21", NULL));
493:   PetscCall(PetscOptionsDeprecated("-tao_cmonitor", "-tao_monitor_constraint_norm", "3.21", NULL));
494:   PetscCall(PetscOptionsDeprecated("-tao_gmonitor", "-tao_monitor_globalization", "3.21", NULL));
495:   PetscCall(PetscOptionsDeprecated("-tao_draw_solution", "-tao_monitor_solution_draw", "3.21", NULL));
496:   PetscCall(PetscOptionsDeprecated("-tao_draw_gradient", "-tao_monitor_gradient_draw", "3.21", NULL));
497:   PetscCall(PetscOptionsDeprecated("-tao_draw_step", "-tao_monitor_step_draw", "3.21", NULL));

499:   PetscCall(PetscOptionsBool("-tao_converged_reason", "Print reason for Tao converged", "TaoSolve", tao->printreason, &tao->printreason, NULL));

501:   PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_solution", "View solution vector after each iteration", "TaoMonitorSolution", TaoMonitorSolution));
502:   PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_gradient", "View gradient vector for each iteration", "TaoMonitorGradient", TaoMonitorGradient));

504:   PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_step", "View step vector after each iteration", "TaoMonitorStep", TaoMonitorStep));
505:   PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_residual", "View least-squares residual vector after each iteration", "TaoMonitorResidual", TaoMonitorResidual));
506:   PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor", "Use the default convergence monitor", "TaoMonitorDefault", TaoMonitorDefault));
507:   PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_globalization", "Use the convergence monitor with extra globalization info", "TaoMonitorGlobalization", TaoMonitorGlobalization));
508:   PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_short", "Use the short convergence monitor", "TaoMonitorDefaultShort", TaoMonitorDefaultShort));
509:   PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_constraint_norm", "Use the default convergence monitor with constraint norm", "TaoMonitorConstraintNorm", TaoMonitorConstraintNorm));

511:   flg = PETSC_FALSE;
512:   PetscCall(PetscOptionsDeprecated("-tao_cancelmonitors", "-tao_monitor_cancel", "3.21", NULL));
513:   PetscCall(PetscOptionsBool("-tao_monitor_cancel", "cancel all monitors and call any registered destroy routines", "TaoMonitorCancel", flg, &flg, NULL));
514:   if (flg) PetscCall(TaoMonitorCancel(tao));

516:   flg = PETSC_FALSE;
517:   PetscCall(PetscOptionsBool("-tao_monitor_solution_draw", "Plot solution vector at each iteration", "TaoMonitorSet", flg, &flg, NULL));
518:   if (flg) {
519:     TaoMonitorDrawCtx drawctx;
520:     PetscInt          howoften = 1;
521:     PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &drawctx));
522:     PetscCall(TaoMonitorSet(tao, TaoMonitorSolutionDraw, drawctx, (PetscCtxDestroyFn *)TaoMonitorDrawCtxDestroy));
523:   }

525:   flg = PETSC_FALSE;
526:   PetscCall(PetscOptionsBool("-tao_monitor_step_draw", "Plots step at each iteration", "TaoMonitorSet", flg, &flg, NULL));
527:   if (flg) PetscCall(TaoMonitorSet(tao, TaoMonitorStepDraw, NULL, NULL));

529:   flg = PETSC_FALSE;
530:   PetscCall(PetscOptionsBool("-tao_monitor_gradient_draw", "plots gradient at each iteration", "TaoMonitorSet", flg, &flg, NULL));
531:   if (flg) {
532:     TaoMonitorDrawCtx drawctx;
533:     PetscInt          howoften = 1;
534:     PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &drawctx));
535:     PetscCall(TaoMonitorSet(tao, TaoMonitorGradientDraw, drawctx, (PetscCtxDestroyFn *)TaoMonitorDrawCtxDestroy));
536:   }
537:   flg = PETSC_FALSE;
538:   PetscCall(PetscOptionsBool("-tao_fd_gradient", "compute gradient using finite differences", "TaoDefaultComputeGradient", flg, &flg, NULL));
539:   if (flg) PetscCall(TaoSetGradient(tao, NULL, TaoDefaultComputeGradient, NULL));
540:   flg = PETSC_FALSE;
541:   PetscCall(PetscOptionsBool("-tao_fd_hessian", "compute Hessian using finite differences", "TaoDefaultComputeHessian", flg, &flg, NULL));
542:   if (flg) {
543:     Mat H;

545:     PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &H));
546:     PetscCall(MatSetType(H, MATAIJ));
547:     PetscCall(TaoSetHessian(tao, H, H, TaoDefaultComputeHessian, NULL));
548:     PetscCall(MatDestroy(&H));
549:   }
550:   flg = PETSC_FALSE;
551:   PetscCall(PetscOptionsBool("-tao_mf_hessian", "compute matrix-free Hessian using finite differences", "TaoDefaultComputeHessianMFFD", flg, &flg, NULL));
552:   if (flg) {
553:     Mat H;

555:     PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &H));
556:     PetscCall(TaoSetHessian(tao, H, H, TaoDefaultComputeHessianMFFD, NULL));
557:     PetscCall(MatDestroy(&H));
558:   }
559:   PetscCall(PetscOptionsBool("-tao_recycle_history", "enable recycling/re-using information from the previous TaoSolve() call for some algorithms", "TaoSetRecycleHistory", flg, &flg, &found));
560:   if (found) PetscCall(TaoSetRecycleHistory(tao, flg));
561:   PetscCall(PetscOptionsEnum("-tao_subset_type", "subset type", "", TaoSubSetTypes, (PetscEnum)tao->subset_type, (PetscEnum *)&tao->subset_type, NULL));

563:   if (tao->ksp) {
564:     PetscCall(PetscOptionsBool("-tao_ksp_ew", "Use Eisentat-Walker linear system convergence test", "TaoKSPSetUseEW", tao->ksp_ewconv, &tao->ksp_ewconv, NULL));
565:     PetscCall(TaoKSPSetUseEW(tao, tao->ksp_ewconv));
566:   }

568:   PetscTryTypeMethod(tao, setfromoptions, PetscOptionsObject);

570:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
571:   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)tao, PetscOptionsObject));
572:   PetscOptionsEnd();

574:   if (tao->linesearch) PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
575:   PetscFunctionReturn(PETSC_SUCCESS);
576: }

578: /*@
579:   TaoViewFromOptions - View a `Tao` object based on values in the options database

581:   Collective

583:   Input Parameters:
584: + A    - the  `Tao` context
585: . obj  - Optional object that provides the prefix for the options database
586: - name - command line option

588:   Level: intermediate

590: .seealso: [](ch_tao), `Tao`, `TaoView`, `PetscObjectViewFromOptions()`, `TaoCreate()`
591: @*/
592: PetscErrorCode TaoViewFromOptions(Tao A, PetscObject obj, const char name[])
593: {
594:   PetscFunctionBegin;
596:   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
597:   PetscFunctionReturn(PETSC_SUCCESS);
598: }

600: /*@
601:   TaoView - Prints information about the `Tao` object

603:   Collective

605:   Input Parameters:
606: + tao    - the `Tao` context
607: - viewer - visualization context

609:   Options Database Key:
610: . -tao_view - Calls `TaoView()` at the end of `TaoSolve()`

612:   Level: beginner

614:   Notes:
615:   The available visualization contexts include
616: +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
617: -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
618:   output where only the first processor opens
619:   the file.  All other processors send their
620:   data to the first processor to print.

622: .seealso: [](ch_tao), `Tao`, `PetscViewerASCIIOpen()`
623: @*/
624: PetscErrorCode TaoView(Tao tao, PetscViewer viewer)
625: {
626:   PetscBool isascii, isstring;
627:   TaoType   type;

629:   PetscFunctionBegin;
631:   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(((PetscObject)tao)->comm, &viewer));
633:   PetscCheckSameComm(tao, 1, viewer, 2);

635:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
636:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
637:   if (isascii) {
638:     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)tao, viewer));

640:     PetscCall(PetscViewerASCIIPushTab(viewer));
641:     PetscTryTypeMethod(tao, view, viewer);
642:     if (tao->linesearch) PetscCall(TaoLineSearchView(tao->linesearch, viewer));
643:     if (tao->ksp) {
644:       PetscCall(KSPView(tao->ksp, viewer));
645:       PetscCall(PetscViewerASCIIPrintf(viewer, "total KSP iterations: %" PetscInt_FMT "\n", tao->ksp_tot_its));
646:     }

648:     if (tao->XL || tao->XU) PetscCall(PetscViewerASCIIPrintf(viewer, "Active Set subset type: %s\n", TaoSubSetTypes[tao->subset_type]));

650:     PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: gatol=%g,", (double)tao->gatol));
651:     PetscCall(PetscViewerASCIIPrintf(viewer, " grtol=%g,", (double)tao->grtol));
652:     PetscCall(PetscViewerASCIIPrintf(viewer, " steptol=%g,", (double)tao->steptol));
653:     PetscCall(PetscViewerASCIIPrintf(viewer, " gttol=%g\n", (double)tao->gttol));
654:     PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Function/Gradient:=%g\n", (double)tao->residual));

656:     if (tao->constrained) {
657:       PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances:"));
658:       PetscCall(PetscViewerASCIIPrintf(viewer, " catol=%g,", (double)tao->catol));
659:       PetscCall(PetscViewerASCIIPrintf(viewer, " crtol=%g\n", (double)tao->crtol));
660:       PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Constraints:=%g\n", (double)tao->cnorm));
661:     }

663:     if (tao->trust < tao->steptol) {
664:       PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: steptol=%g\n", (double)tao->steptol));
665:       PetscCall(PetscViewerASCIIPrintf(viewer, "Final trust region radius:=%g\n", (double)tao->trust));
666:     }

668:     if (tao->fmin > -1.e25) PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: function minimum=%g\n", (double)tao->fmin));
669:     PetscCall(PetscViewerASCIIPrintf(viewer, "Objective value=%g\n", (double)tao->fc));

671:     PetscCall(PetscViewerASCIIPrintf(viewer, "total number of iterations=%" PetscInt_FMT ",          ", tao->niter));
672:     PetscCall(PetscViewerASCIIPrintf(viewer, "              (max: %" PetscInt_FMT ")\n", tao->max_it));

674:     if (tao->nfuncs > 0) {
675:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function evaluations=%" PetscInt_FMT ",", tao->nfuncs));
676:       if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: unlimited)\n"));
677:       else PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: %" PetscInt_FMT ")\n", tao->max_funcs));
678:     }
679:     if (tao->ngrads > 0) {
680:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of gradient evaluations=%" PetscInt_FMT ",", tao->ngrads));
681:       if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: unlimited)\n"));
682:       else PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: %" PetscInt_FMT ")\n", tao->max_funcs));
683:     }
684:     if (tao->nfuncgrads > 0) {
685:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function/gradient evaluations=%" PetscInt_FMT ",", tao->nfuncgrads));
686:       if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, "    (max: unlimited)\n"));
687:       else PetscCall(PetscViewerASCIIPrintf(viewer, "    (max: %" PetscInt_FMT ")\n", tao->max_funcs));
688:     }
689:     if (tao->nhess > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Hessian evaluations=%" PetscInt_FMT "\n", tao->nhess));
690:     if (tao->nconstraints > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of constraint function evaluations=%" PetscInt_FMT "\n", tao->nconstraints));
691:     if (tao->njac > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Jacobian evaluations=%" PetscInt_FMT "\n", tao->njac));

693:     if (tao->reason > 0) {
694:       PetscCall(PetscViewerASCIIPrintf(viewer, "Solution converged: "));
695:       switch (tao->reason) {
696:       case TAO_CONVERGED_GATOL:
697:         PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)|| <= gatol\n"));
698:         break;
699:       case TAO_CONVERGED_GRTOL:
700:         PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)||/|f(X)| <= grtol\n"));
701:         break;
702:       case TAO_CONVERGED_GTTOL:
703:         PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)||/||g(X0)|| <= gttol\n"));
704:         break;
705:       case TAO_CONVERGED_STEPTOL:
706:         PetscCall(PetscViewerASCIIPrintf(viewer, " Steptol -- step size small\n"));
707:         break;
708:       case TAO_CONVERGED_MINF:
709:         PetscCall(PetscViewerASCIIPrintf(viewer, " Minf --  f < fmin\n"));
710:         break;
711:       case TAO_CONVERGED_USER:
712:         PetscCall(PetscViewerASCIIPrintf(viewer, " User Terminated\n"));
713:         break;
714:       default:
715:         PetscCall(PetscViewerASCIIPrintf(viewer, " %d\n", tao->reason));
716:         break;
717:       }
718:     } else if (tao->reason == TAO_CONTINUE_ITERATING) {
719:       PetscCall(PetscViewerASCIIPrintf(viewer, "Solver never run\n"));
720:     } else {
721:       PetscCall(PetscViewerASCIIPrintf(viewer, "Solver failed: "));
722:       switch (tao->reason) {
723:       case TAO_DIVERGED_MAXITS:
724:         PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum Iterations\n"));
725:         break;
726:       case TAO_DIVERGED_NAN:
727:         PetscCall(PetscViewerASCIIPrintf(viewer, " NAN or infinity encountered\n"));
728:         break;
729:       case TAO_DIVERGED_MAXFCN:
730:         PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum Function Evaluations\n"));
731:         break;
732:       case TAO_DIVERGED_LS_FAILURE:
733:         PetscCall(PetscViewerASCIIPrintf(viewer, " Line Search Failure\n"));
734:         break;
735:       case TAO_DIVERGED_TR_REDUCTION:
736:         PetscCall(PetscViewerASCIIPrintf(viewer, " Trust Region too small\n"));
737:         break;
738:       case TAO_DIVERGED_USER:
739:         PetscCall(PetscViewerASCIIPrintf(viewer, " User Terminated\n"));
740:         break;
741:       default:
742:         PetscCall(PetscViewerASCIIPrintf(viewer, " %d\n", tao->reason));
743:         break;
744:       }
745:     }
746:     PetscCall(PetscViewerASCIIPopTab(viewer));
747:   } else if (isstring) {
748:     PetscCall(TaoGetType(tao, &type));
749:     PetscCall(PetscViewerStringSPrintf(viewer, " %-3.3s", type));
750:   }
751:   PetscFunctionReturn(PETSC_SUCCESS);
752: }

754: /*@
755:   TaoSetRecycleHistory - Sets the boolean flag to enable/disable re-using
756:   iterate information from the previous `TaoSolve()`. This feature is disabled by
757:   default.

759:   Logically Collective

761:   Input Parameters:
762: + tao     - the `Tao` context
763: - recycle - boolean flag

765:   Options Database Key:
766: . -tao_recycle_history <true,false> - reuse the history

768:   Level: intermediate

770:   Notes:
771:   For conjugate gradient methods (`TAOBNCG`), this re-uses the latest search direction
772:   from the previous `TaoSolve()` call when computing the first search direction in a
773:   new solution. By default, CG methods set the first search direction to the
774:   negative gradient.

776:   For quasi-Newton family of methods (`TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`), this re-uses
777:   the accumulated quasi-Newton Hessian approximation from the previous `TaoSolve()`
778:   call. By default, QN family of methods reset the initial Hessian approximation to
779:   the identity matrix.

781:   For any other algorithm, this setting has no effect.

783: .seealso: [](ch_tao), `Tao`, `TaoGetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
784: @*/
785: PetscErrorCode TaoSetRecycleHistory(Tao tao, PetscBool recycle)
786: {
787:   PetscFunctionBegin;
790:   tao->recycle = recycle;
791:   PetscFunctionReturn(PETSC_SUCCESS);
792: }

794: /*@
795:   TaoGetRecycleHistory - Retrieve the boolean flag for re-using iterate information
796:   from the previous `TaoSolve()`. This feature is disabled by default.

798:   Logically Collective

800:   Input Parameter:
801: . tao - the `Tao` context

803:   Output Parameter:
804: . recycle - boolean flag

806:   Level: intermediate

808: .seealso: [](ch_tao), `Tao`, `TaoSetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
809: @*/
810: PetscErrorCode TaoGetRecycleHistory(Tao tao, PetscBool *recycle)
811: {
812:   PetscFunctionBegin;
814:   PetscAssertPointer(recycle, 2);
815:   *recycle = tao->recycle;
816:   PetscFunctionReturn(PETSC_SUCCESS);
817: }

819: /*@
820:   TaoSetTolerances - Sets parameters used in `TaoSolve()` convergence tests

822:   Logically Collective

824:   Input Parameters:
825: + tao   - the `Tao` context
826: . gatol - stop if norm of gradient is less than this
827: . grtol - stop if relative norm of gradient is less than this
828: - gttol - stop if norm of gradient is reduced by this factor

830:   Options Database Keys:
831: + -tao_gatol <gatol> - Sets gatol
832: . -tao_grtol <grtol> - Sets grtol
833: - -tao_gttol <gttol> - Sets gttol

835:   Stopping Criteria\:
836: .vb
837:   ||g(X)||                            <= gatol
838:   ||g(X)|| / |f(X)|                   <= grtol
839:   ||g(X)|| / ||g(X0)||                <= gttol
840: .ve

842:   Level: beginner

844:   Notes:
845:   Use `PETSC_CURRENT` to leave one or more tolerances unchanged.

847:   Use `PETSC_DETERMINE` to set one or more tolerances to their values when the `tao`object's type was set

849:   Fortran Note:
850:   Use `PETSC_CURRENT_REAL` or `PETSC_DETERMINE_REAL`

852: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`
853: @*/
854: PetscErrorCode TaoSetTolerances(Tao tao, PetscReal gatol, PetscReal grtol, PetscReal gttol)
855: {
856:   PetscFunctionBegin;

862:   if (gatol == (PetscReal)PETSC_DETERMINE) {
863:     tao->gatol = tao->default_gatol;
864:   } else if (gatol != (PetscReal)PETSC_CURRENT) {
865:     PetscCheck(gatol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative gatol not allowed");
866:     tao->gatol = gatol;
867:   }

869:   if (grtol == (PetscReal)PETSC_DETERMINE) {
870:     tao->grtol = tao->default_grtol;
871:   } else if (grtol != (PetscReal)PETSC_CURRENT) {
872:     PetscCheck(grtol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative grtol not allowed");
873:     tao->grtol = grtol;
874:   }

876:   if (gttol == (PetscReal)PETSC_DETERMINE) {
877:     tao->gttol = tao->default_gttol;
878:   } else if (gttol != (PetscReal)PETSC_CURRENT) {
879:     PetscCheck(gttol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative gttol not allowed");
880:     tao->gttol = gttol;
881:   }
882:   PetscFunctionReturn(PETSC_SUCCESS);
883: }

885: /*@
886:   TaoSetConstraintTolerances - Sets constraint tolerance parameters used in `TaoSolve()` convergence tests

888:   Logically Collective

890:   Input Parameters:
891: + tao   - the `Tao` context
892: . catol - absolute constraint tolerance, constraint norm must be less than `catol` for used for `gatol` convergence criteria
893: - crtol - relative constraint tolerance, constraint norm must be less than `crtol` for used for `gatol`, `gttol` convergence criteria

895:   Options Database Keys:
896: + -tao_catol <catol> - Sets catol
897: - -tao_crtol <crtol> - Sets crtol

899:   Level: intermediate

901:   Notes:
902:   Use `PETSC_CURRENT` to leave one or tolerance unchanged.

904:   Use `PETSC_DETERMINE` to set one or more tolerances to their values when the `tao` object's type was set

906:   Fortran Note:
907:   Use `PETSC_CURRENT_REAL` or `PETSC_DETERMINE_REAL`

909: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`, `TaoGetConstraintTolerances()`, `TaoSetTolerances()`
910: @*/
911: PetscErrorCode TaoSetConstraintTolerances(Tao tao, PetscReal catol, PetscReal crtol)
912: {
913:   PetscFunctionBegin;

918:   if (catol == (PetscReal)PETSC_DETERMINE) {
919:     tao->catol = tao->default_catol;
920:   } else if (catol != (PetscReal)PETSC_CURRENT) {
921:     PetscCheck(catol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative catol not allowed");
922:     tao->catol = catol;
923:   }

925:   if (crtol == (PetscReal)PETSC_DETERMINE) {
926:     tao->crtol = tao->default_crtol;
927:   } else if (crtol != (PetscReal)PETSC_CURRENT) {
928:     PetscCheck(crtol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative crtol not allowed");
929:     tao->crtol = crtol;
930:   }
931:   PetscFunctionReturn(PETSC_SUCCESS);
932: }

934: /*@
935:   TaoGetConstraintTolerances - Gets constraint tolerance parameters used in `TaoSolve()` convergence tests

937:   Not Collective

939:   Input Parameter:
940: . tao - the `Tao` context

942:   Output Parameters:
943: + catol - absolute constraint tolerance, constraint norm must be less than `catol` for used for `gatol` convergence criteria
944: - crtol - relative constraint tolerance, constraint norm must be less than `crtol` for used for `gatol`, `gttol` convergence criteria

946:   Level: intermediate

948: .seealso: [](ch_tao), `Tao`, `TaoConvergedReasons`,`TaoGetTolerances()`, `TaoSetTolerances()`, `TaoSetConstraintTolerances()`
949: @*/
950: PetscErrorCode TaoGetConstraintTolerances(Tao tao, PetscReal *catol, PetscReal *crtol)
951: {
952:   PetscFunctionBegin;
954:   if (catol) *catol = tao->catol;
955:   if (crtol) *crtol = tao->crtol;
956:   PetscFunctionReturn(PETSC_SUCCESS);
957: }

959: /*@
960:   TaoSetFunctionLowerBound - Sets a 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.

964:   Logically Collective

966:   Input Parameters:
967: + tao  - the Tao solver context
968: - fmin - the tolerance

970:   Options Database Key:
971: . -tao_fmin <fmin> - sets the minimum function value

973:   Level: intermediate

975: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetTolerances()`
976: @*/
977: PetscErrorCode TaoSetFunctionLowerBound(Tao tao, PetscReal fmin)
978: {
979:   PetscFunctionBegin;
982:   tao->fmin = fmin;
983:   PetscFunctionReturn(PETSC_SUCCESS);
984: }

986: /*@
987:   TaoGetFunctionLowerBound - Gets the bound on the solution objective value.
988:   When an approximate solution with an objective value below this number
989:   has been found, the solver will terminate.

991:   Not Collective

993:   Input Parameter:
994: . tao - the `Tao` solver context

996:   Output Parameter:
997: . fmin - the minimum function value

999:   Level: intermediate

1001: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetFunctionLowerBound()`
1002: @*/
1003: PetscErrorCode TaoGetFunctionLowerBound(Tao tao, PetscReal *fmin)
1004: {
1005:   PetscFunctionBegin;
1007:   PetscAssertPointer(fmin, 2);
1008:   *fmin = tao->fmin;
1009:   PetscFunctionReturn(PETSC_SUCCESS);
1010: }

1012: /*@
1013:   TaoSetMaximumFunctionEvaluations - Sets a maximum number of function evaluations allowed for a `TaoSolve()`.

1015:   Logically Collective

1017:   Input Parameters:
1018: + tao  - the `Tao` solver context
1019: - nfcn - the maximum number of function evaluations (>=0), use `PETSC_UNLIMITED` to have no bound

1021:   Options Database Key:
1022: . -tao_max_funcs <nfcn> - sets the maximum number of function evaluations

1024:   Level: intermediate

1026:   Note:
1027:   Use `PETSC_DETERMINE` to use the default maximum number of function evaluations that was set when the object type was set.

1029:   Developer Note:
1030:   Deprecated support for an unlimited number of function evaluations by passing a negative value.

1032: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoSetMaximumIterations()`
1033: @*/
1034: PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao tao, PetscInt nfcn)
1035: {
1036:   PetscFunctionBegin;
1039:   if (nfcn == PETSC_DETERMINE) {
1040:     tao->max_funcs = tao->default_max_funcs;
1041:   } else if (nfcn == PETSC_UNLIMITED || nfcn < 0) {
1042:     tao->max_funcs = PETSC_UNLIMITED;
1043:   } else {
1044:     PetscCheck(nfcn >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of function evaluations  must be positive");
1045:     tao->max_funcs = nfcn;
1046:   }
1047:   PetscFunctionReturn(PETSC_SUCCESS);
1048: }

1050: /*@
1051:   TaoGetMaximumFunctionEvaluations - Gets a maximum number of function evaluations allowed for a `TaoSolve()`

1053:   Logically Collective

1055:   Input Parameter:
1056: . tao - the `Tao` solver context

1058:   Output Parameter:
1059: . nfcn - the maximum number of function evaluations

1061:   Level: intermediate

1063: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1064: @*/
1065: PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao tao, PetscInt *nfcn)
1066: {
1067:   PetscFunctionBegin;
1069:   PetscAssertPointer(nfcn, 2);
1070:   *nfcn = tao->max_funcs;
1071:   PetscFunctionReturn(PETSC_SUCCESS);
1072: }

1074: /*@
1075:   TaoGetCurrentFunctionEvaluations - Get current number of function evaluations used by a `Tao` object

1077:   Not Collective

1079:   Input Parameter:
1080: . tao - the `Tao` solver context

1082:   Output Parameter:
1083: . nfuncs - the current number of function evaluations (maximum between gradient and function evaluations)

1085:   Level: intermediate

1087: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1088: @*/
1089: PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao tao, PetscInt *nfuncs)
1090: {
1091:   PetscFunctionBegin;
1093:   PetscAssertPointer(nfuncs, 2);
1094:   *nfuncs = PetscMax(tao->nfuncs, tao->nfuncgrads);
1095:   PetscFunctionReturn(PETSC_SUCCESS);
1096: }

1098: /*@
1099:   TaoSetMaximumIterations - Sets a maximum number of iterates to be used in `TaoSolve()`

1101:   Logically Collective

1103:   Input Parameters:
1104: + tao    - the `Tao` solver context
1105: - maxits - the maximum number of iterates (>=0), use `PETSC_UNLIMITED` to have no bound

1107:   Options Database Key:
1108: . -tao_max_it <its> - sets the maximum number of iterations

1110:   Level: intermediate

1112:   Note:
1113:   Use `PETSC_DETERMINE` to use the default maximum number of iterations that was set when the object's type was set.

1115:   Developer Note:
1116:   DeprAlso accepts the deprecated negative values to indicate no limit

1118: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoSetMaximumFunctionEvaluations()`
1119: @*/
1120: PetscErrorCode TaoSetMaximumIterations(Tao tao, PetscInt maxits)
1121: {
1122:   PetscFunctionBegin;
1125:   if (maxits == PETSC_DETERMINE) {
1126:     tao->max_it = tao->default_max_it;
1127:   } else if (maxits == PETSC_UNLIMITED) {
1128:     tao->max_it = PETSC_INT_MAX;
1129:   } else {
1130:     PetscCheck(maxits > 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations must be positive");
1131:     tao->max_it = maxits;
1132:   }
1133:   PetscFunctionReturn(PETSC_SUCCESS);
1134: }

1136: /*@
1137:   TaoGetMaximumIterations - Gets a maximum number of iterates that will be used

1139:   Not Collective

1141:   Input Parameter:
1142: . tao - the `Tao` solver context

1144:   Output Parameter:
1145: . maxits - the maximum number of iterates

1147:   Level: intermediate

1149: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumIterations()`, `TaoGetMaximumFunctionEvaluations()`
1150: @*/
1151: PetscErrorCode TaoGetMaximumIterations(Tao tao, PetscInt *maxits)
1152: {
1153:   PetscFunctionBegin;
1155:   PetscAssertPointer(maxits, 2);
1156:   *maxits = tao->max_it;
1157:   PetscFunctionReturn(PETSC_SUCCESS);
1158: }

1160: /*@
1161:   TaoSetInitialTrustRegionRadius - Sets the initial trust region radius.

1163:   Logically Collective

1165:   Input Parameters:
1166: + tao    - a `Tao` optimization solver
1167: - radius - the trust region radius

1169:   Options Database Key:
1170: . -tao_trust0 <t0> - sets initial trust region radius

1172:   Level: intermediate

1174:   Note:
1175:   Use `PETSC_DETERMINE` to use the default radius that was set when the object's type was set.

1177: .seealso: [](ch_tao), `Tao`, `TaoGetTrustRegionRadius()`, `TaoSetTrustRegionTolerance()`, `TAONTR`
1178: @*/
1179: PetscErrorCode TaoSetInitialTrustRegionRadius(Tao tao, PetscReal radius)
1180: {
1181:   PetscFunctionBegin;
1184:   if (radius == PETSC_DETERMINE) {
1185:     tao->trust0 = tao->default_trust0;
1186:   } else {
1187:     PetscCheck(radius > 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Radius must be positive");
1188:     tao->trust0 = radius;
1189:   }
1190:   PetscFunctionReturn(PETSC_SUCCESS);
1191: }

1193: /*@
1194:   TaoGetInitialTrustRegionRadius - Gets the initial trust region radius.

1196:   Not Collective

1198:   Input Parameter:
1199: . tao - a `Tao` optimization solver

1201:   Output Parameter:
1202: . radius - the trust region radius

1204:   Level: intermediate

1206: .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetCurrentTrustRegionRadius()`, `TAONTR`
1207: @*/
1208: PetscErrorCode TaoGetInitialTrustRegionRadius(Tao tao, PetscReal *radius)
1209: {
1210:   PetscFunctionBegin;
1212:   PetscAssertPointer(radius, 2);
1213:   *radius = tao->trust0;
1214:   PetscFunctionReturn(PETSC_SUCCESS);
1215: }

1217: /*@
1218:   TaoGetCurrentTrustRegionRadius - Gets the current trust region radius.

1220:   Not Collective

1222:   Input Parameter:
1223: . tao - a `Tao` optimization solver

1225:   Output Parameter:
1226: . radius - the trust region radius

1228:   Level: intermediate

1230: .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetInitialTrustRegionRadius()`, `TAONTR`
1231: @*/
1232: PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao tao, PetscReal *radius)
1233: {
1234:   PetscFunctionBegin;
1236:   PetscAssertPointer(radius, 2);
1237:   *radius = tao->trust;
1238:   PetscFunctionReturn(PETSC_SUCCESS);
1239: }

1241: /*@
1242:   TaoGetTolerances - gets the current values of some tolerances used for the convergence testing of `TaoSolve()`

1244:   Not Collective

1246:   Input Parameter:
1247: . tao - the `Tao` context

1249:   Output Parameters:
1250: + gatol - stop if norm of gradient is less than this
1251: . grtol - stop if relative norm of gradient is less than this
1252: - gttol - stop if norm of gradient is reduced by a this factor

1254:   Level: intermediate

1256:   Note:
1257:   `NULL` can be used as an argument if not all tolerances values are needed

1259: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`
1260: @*/
1261: PetscErrorCode TaoGetTolerances(Tao tao, PetscReal *gatol, PetscReal *grtol, PetscReal *gttol)
1262: {
1263:   PetscFunctionBegin;
1265:   if (gatol) *gatol = tao->gatol;
1266:   if (grtol) *grtol = tao->grtol;
1267:   if (gttol) *gttol = tao->gttol;
1268:   PetscFunctionReturn(PETSC_SUCCESS);
1269: }

1271: /*@
1272:   TaoGetKSP - Gets the linear solver used by the optimization solver.

1274:   Not Collective

1276:   Input Parameter:
1277: . tao - the `Tao` solver

1279:   Output Parameter:
1280: . ksp - the `KSP` linear solver used in the optimization solver

1282:   Level: intermediate

1284: .seealso: [](ch_tao), `Tao`, `KSP`
1285: @*/
1286: PetscErrorCode TaoGetKSP(Tao tao, KSP *ksp)
1287: {
1288:   PetscFunctionBegin;
1290:   PetscAssertPointer(ksp, 2);
1291:   *ksp = tao->ksp;
1292:   PetscFunctionReturn(PETSC_SUCCESS);
1293: }

1295: /*@
1296:   TaoGetLinearSolveIterations - Gets the total number of linear iterations
1297:   used by the `Tao` solver

1299:   Not Collective

1301:   Input Parameter:
1302: . tao - the `Tao` context

1304:   Output Parameter:
1305: . lits - number of linear iterations

1307:   Level: intermediate

1309:   Note:
1310:   This counter is reset to zero for each successive call to `TaoSolve()`

1312: .seealso: [](ch_tao), `Tao`, `TaoGetKSP()`
1313: @*/
1314: PetscErrorCode TaoGetLinearSolveIterations(Tao tao, PetscInt *lits)
1315: {
1316:   PetscFunctionBegin;
1318:   PetscAssertPointer(lits, 2);
1319:   *lits = tao->ksp_tot_its;
1320:   PetscFunctionReturn(PETSC_SUCCESS);
1321: }

1323: /*@
1324:   TaoGetLineSearch - Gets the line search used by the optimization solver.

1326:   Not Collective

1328:   Input Parameter:
1329: . tao - the `Tao` solver

1331:   Output Parameter:
1332: . ls - the line search used in the optimization solver

1334:   Level: intermediate

1336: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`
1337: @*/
1338: PetscErrorCode TaoGetLineSearch(Tao tao, TaoLineSearch *ls)
1339: {
1340:   PetscFunctionBegin;
1342:   PetscAssertPointer(ls, 2);
1343:   *ls = tao->linesearch;
1344:   PetscFunctionReturn(PETSC_SUCCESS);
1345: }

1347: /*@
1348:   TaoAddLineSearchCounts - Adds the number of function evaluations spent
1349:   in the line search to the running total.

1351:   Input Parameters:
1352: . tao - the `Tao` solver

1354:   Level: developer

1356: .seealso: [](ch_tao), `Tao`, `TaoGetLineSearch()`, `TaoLineSearchApply()`
1357: @*/
1358: PetscErrorCode TaoAddLineSearchCounts(Tao tao)
1359: {
1360:   PetscBool flg;
1361:   PetscInt  nfeval, ngeval, nfgeval;

1363:   PetscFunctionBegin;
1365:   if (tao->linesearch) {
1366:     PetscCall(TaoLineSearchIsUsingTaoRoutines(tao->linesearch, &flg));
1367:     if (!flg) {
1368:       PetscCall(TaoLineSearchGetNumberFunctionEvaluations(tao->linesearch, &nfeval, &ngeval, &nfgeval));
1369:       tao->nfuncs += nfeval;
1370:       tao->ngrads += ngeval;
1371:       tao->nfuncgrads += nfgeval;
1372:     }
1373:   }
1374:   PetscFunctionReturn(PETSC_SUCCESS);
1375: }

1377: /*@
1378:   TaoGetSolution - Returns the vector with the current solution from the `Tao` object

1380:   Not Collective

1382:   Input Parameter:
1383: . tao - the `Tao` context

1385:   Output Parameter:
1386: . X - the current solution

1388:   Level: intermediate

1390:   Note:
1391:   The returned vector will be the same object that was passed into `TaoSetSolution()`

1393: .seealso: [](ch_tao), `Tao`, `TaoSetSolution()`, `TaoSolve()`
1394: @*/
1395: PetscErrorCode TaoGetSolution(Tao tao, Vec *X)
1396: {
1397:   PetscFunctionBegin;
1399:   PetscAssertPointer(X, 2);
1400:   *X = tao->solution;
1401:   PetscFunctionReturn(PETSC_SUCCESS);
1402: }

1404: /*@
1405:   TaoResetStatistics - Initialize the statistics collected by the `Tao` object.
1406:   These statistics include the iteration number, residual norms, and convergence status.
1407:   This routine gets called before solving each optimization problem.

1409:   Collective

1411:   Input Parameter:
1412: . tao - the `Tao` context

1414:   Level: developer

1416: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
1417: @*/
1418: PetscErrorCode TaoResetStatistics(Tao tao)
1419: {
1420:   PetscFunctionBegin;
1422:   tao->niter        = 0;
1423:   tao->nfuncs       = 0;
1424:   tao->nfuncgrads   = 0;
1425:   tao->ngrads       = 0;
1426:   tao->nhess        = 0;
1427:   tao->njac         = 0;
1428:   tao->nconstraints = 0;
1429:   tao->ksp_its      = 0;
1430:   tao->ksp_tot_its  = 0;
1431:   tao->reason       = TAO_CONTINUE_ITERATING;
1432:   tao->residual     = 0.0;
1433:   tao->cnorm        = 0.0;
1434:   tao->step         = 0.0;
1435:   tao->lsflag       = PETSC_FALSE;
1436:   if (tao->hist_reset) tao->hist_len = 0;
1437:   PetscFunctionReturn(PETSC_SUCCESS);
1438: }

1440: /*@C
1441:   TaoSetUpdate - Sets the general-purpose update function called
1442:   at the beginning of every iteration of the optimization algorithm. Called after the new solution and the gradient
1443:   is determined, but before the Hessian is computed (if applicable).

1445:   Logically Collective

1447:   Input Parameters:
1448: + tao  - The `Tao` solver
1449: . func - The function
1450: - ctx  - The update function context

1452:   Calling sequence of `func`:
1453: + tao - The optimizer context
1454: . it  - The current iteration index
1455: - ctx - The update context

1457:   Level: advanced

1459:   Notes:
1460:   Users can modify the gradient direction or any other vector associated to the specific solver used.
1461:   The objective function value is always recomputed after a call to the update hook.

1463: .seealso: [](ch_tao), `Tao`, `TaoSolve()`
1464: @*/
1465: PetscErrorCode TaoSetUpdate(Tao tao, PetscErrorCode (*func)(Tao tao, PetscInt it, PetscCtx ctx), PetscCtx ctx)
1466: {
1467:   PetscFunctionBegin;
1469:   tao->ops->update = func;
1470:   tao->user_update = ctx;
1471:   PetscFunctionReturn(PETSC_SUCCESS);
1472: }

1474: /*@C
1475:   TaoSetConvergenceTest - Sets the function that is to be used to test
1476:   for convergence of the iterative minimization solution.  The new convergence
1477:   testing routine will replace Tao's default convergence test.

1479:   Logically Collective

1481:   Input Parameters:
1482: + tao  - the `Tao` object
1483: . conv - the routine to test for convergence
1484: - ctx  - [optional] context for private data for the convergence routine (may be `NULL`)

1486:   Calling sequence of `conv`:
1487: + tao - the `Tao` object
1488: - ctx - [optional] convergence context

1490:   Level: advanced

1492:   Note:
1493:   The new convergence testing routine should call `TaoSetConvergedReason()`.

1495: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetConvergedReason()`, `TaoGetSolutionStatus()`, `TaoGetTolerances()`, `TaoMonitorSet()`
1496: @*/
1497: PetscErrorCode TaoSetConvergenceTest(Tao tao, PetscErrorCode (*conv)(Tao tao, PetscCtx ctx), PetscCtx ctx)
1498: {
1499:   PetscFunctionBegin;
1501:   tao->ops->convergencetest = conv;
1502:   tao->cnvP                 = ctx;
1503:   PetscFunctionReturn(PETSC_SUCCESS);
1504: }

1506: /*@C
1507:   TaoMonitorSet - Sets an additional function that is to be used at every
1508:   iteration of the solver to display the iteration's
1509:   progress.

1511:   Logically Collective

1513:   Input Parameters:
1514: + tao  - the `Tao` solver context
1515: . func - monitoring routine
1516: . ctx  - [optional] user-defined context for private data for the monitor routine (may be `NULL`)
1517: - dest - [optional] function to destroy the context when the `Tao` is destroyed, see `PetscCtxDestroyFn` for the calling sequence

1519:   Calling sequence of `func`:
1520: + tao - the `Tao` solver context
1521: - ctx - [optional] monitoring context

1523:   Level: intermediate

1525:   Notes:
1526:   See `TaoSetFromOptions()` for a monitoring options.

1528:   Several different monitoring routines may be set by calling
1529:   `TaoMonitorSet()` multiple times; all will be called in the
1530:   order in which they were set.

1532:   Fortran Notes:
1533:   Only one monitor function may be set

1535: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoMonitorDefault()`, `TaoMonitorCancel()`, `TaoView()`, `PetscCtxDestroyFn`
1536: @*/
1537: PetscErrorCode TaoMonitorSet(Tao tao, PetscErrorCode (*func)(Tao tao, PetscCtx ctx), PetscCtx ctx, PetscCtxDestroyFn *dest)
1538: {
1539:   PetscFunctionBegin;
1541:   PetscCheck(tao->numbermonitors < MAXTAOMONITORS, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Cannot attach another monitor -- max=%d", MAXTAOMONITORS);
1542:   for (PetscInt i = 0; i < tao->numbermonitors; i++) {
1543:     PetscBool identical;

1545:     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))(PetscVoidFn *)func, ctx, dest, (PetscErrorCode (*)(void))(PetscVoidFn *)tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i], &identical));
1546:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
1547:   }
1548:   tao->monitor[tao->numbermonitors]        = func;
1549:   tao->monitorcontext[tao->numbermonitors] = ctx;
1550:   tao->monitordestroy[tao->numbermonitors] = dest;
1551:   ++tao->numbermonitors;
1552:   PetscFunctionReturn(PETSC_SUCCESS);
1553: }

1555: /*@
1556:   TaoMonitorCancel - Clears all the monitor functions for a `Tao` object.

1558:   Logically Collective

1560:   Input Parameter:
1561: . tao - the `Tao` solver context

1563:   Options Database Key:
1564: . -tao_monitor_cancel - cancels all monitors that have been hardwired
1565:     into a code by calls to `TaoMonitorSet()`, but does not cancel those
1566:     set via the options database

1568:   Level: advanced

1570:   Note:
1571:   There is no way to clear one specific monitor from a `Tao` object.

1573: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1574: @*/
1575: PetscErrorCode TaoMonitorCancel(Tao tao)
1576: {
1577:   PetscInt i;

1579:   PetscFunctionBegin;
1581:   for (i = 0; i < tao->numbermonitors; i++) {
1582:     if (tao->monitordestroy[i]) PetscCall((*tao->monitordestroy[i])(&tao->monitorcontext[i]));
1583:   }
1584:   tao->numbermonitors = 0;
1585:   PetscFunctionReturn(PETSC_SUCCESS);
1586: }

1588: /*@
1589:   TaoMonitorDefault - Default routine for monitoring progress of `TaoSolve()`

1591:   Collective

1593:   Input Parameters:
1594: + tao - the `Tao` context
1595: - vf  - `PetscViewerAndFormat` context

1597:   Options Database Key:
1598: . -tao_monitor - turn on default monitoring

1600:   Level: advanced

1602:   Note:
1603:   This monitor prints the function value and gradient
1604:   norm at each iteration.

1606: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1607: @*/
1608: PetscErrorCode TaoMonitorDefault(Tao tao, PetscViewerAndFormat *vf)
1609: {
1610:   PetscViewer viewer = vf->viewer;
1611:   PetscBool   isascii;
1612:   PetscInt    tabs;

1614:   PetscFunctionBegin;
1616:   if (vf->view_interval > 0 && tao->niter % vf->view_interval) PetscFunctionReturn(PETSC_SUCCESS);

1618:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1619:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1620:   if (isascii) {
1621:     PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));

1623:     PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1624:     if (tao->niter == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
1625:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Iteration information for %s solve.\n", ((PetscObject)tao)->prefix));
1626:       tao->header_printed = PETSC_TRUE;
1627:     }
1628:     PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " TAO,", tao->niter));
1629:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Function value: %g,", (double)tao->fc));
1630:     if (tao->residual >= PETSC_INFINITY) {
1631:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: infinity \n"));
1632:     } else {
1633:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g \n", (double)tao->residual));
1634:     }
1635:     PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1636:   }
1637:   PetscCall(PetscViewerPopFormat(viewer));
1638:   PetscFunctionReturn(PETSC_SUCCESS);
1639: }

1641: /*@
1642:   TaoMonitorGlobalization - Default routine for monitoring progress of `TaoSolve()` with extra detail on the globalization method.

1644:   Collective

1646:   Input Parameters:
1647: + tao - the `Tao` context
1648: - vf  - `PetscViewerAndFormat` context

1650:   Options Database Key:
1651: . -tao_monitor_globalization - turn on monitoring with globalization information

1653:   Level: advanced

1655:   Note:
1656:   This monitor prints the function value and gradient norm at each
1657:   iteration, as well as the step size and trust radius. Note that the
1658:   step size and trust radius may be the same for some algorithms.

1660: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1661: @*/
1662: PetscErrorCode TaoMonitorGlobalization(Tao tao, PetscViewerAndFormat *vf)
1663: {
1664:   PetscViewer viewer = vf->viewer;
1665:   PetscBool   isascii;
1666:   PetscInt    tabs;

1668:   PetscFunctionBegin;
1670:   if (vf->view_interval > 0 && tao->niter % vf->view_interval) PetscFunctionReturn(PETSC_SUCCESS);

1672:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1673:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1674:   if (isascii) {
1675:     PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1676:     PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1677:     if (tao->niter == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
1678:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Iteration information for %s solve.\n", ((PetscObject)tao)->prefix));
1679:       tao->header_printed = PETSC_TRUE;
1680:     }
1681:     PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " TAO,", tao->niter));
1682:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Function value: %g,", (double)tao->fc));
1683:     if (tao->residual >= PETSC_INFINITY) {
1684:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: Inf,"));
1685:     } else {
1686:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g,", (double)tao->residual));
1687:     }
1688:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Step: %g,  Trust: %g\n", (double)tao->step, (double)tao->trust));
1689:     PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1690:   }
1691:   PetscCall(PetscViewerPopFormat(viewer));
1692:   PetscFunctionReturn(PETSC_SUCCESS);
1693: }

1695: /*@
1696:   TaoMonitorDefaultShort - Routine for monitoring progress of `TaoSolve()` that displays fewer digits than `TaoMonitorDefault()`

1698:   Collective

1700:   Input Parameters:
1701: + tao - the `Tao` context
1702: - vf  - `PetscViewerAndFormat` context

1704:   Options Database Key:
1705: . -tao_monitor_short - turn on default short monitoring

1707:   Level: advanced

1709:   Note:
1710:   Same as `TaoMonitorDefault()` except
1711:   it prints fewer digits of the residual as the residual gets smaller.
1712:   This is because the later digits are meaningless and are often
1713:   different on different machines; by using this routine different
1714:   machines will usually generate the same output.

1716: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1717: @*/
1718: PetscErrorCode TaoMonitorDefaultShort(Tao tao, PetscViewerAndFormat *vf)
1719: {
1720:   PetscViewer viewer = vf->viewer;
1721:   PetscBool   isascii;
1722:   PetscInt    tabs;
1723:   PetscReal   gnorm;

1725:   PetscFunctionBegin;
1727:   if (vf->view_interval > 0 && tao->niter % vf->view_interval) PetscFunctionReturn(PETSC_SUCCESS);

1729:   gnorm = tao->residual;
1730:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1731:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1732:   if (isascii) {
1733:     PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1734:     PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1735:     PetscCall(PetscViewerASCIIPrintf(viewer, "iter = %3" PetscInt_FMT ",", tao->niter));
1736:     PetscCall(PetscViewerASCIIPrintf(viewer, " Function value %g,", (double)tao->fc));
1737:     if (gnorm >= PETSC_INFINITY) {
1738:       PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: infinity \n"));
1739:     } else if (gnorm > 1.e-6) {
1740:       PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: %g \n", (double)gnorm));
1741:     } else if (gnorm > 1.e-11) {
1742:       PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: < 1.0e-6 \n"));
1743:     } else {
1744:       PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: < 1.0e-11 \n"));
1745:     }
1746:     PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1747:   }
1748:   PetscCall(PetscViewerPopFormat(viewer));
1749:   PetscFunctionReturn(PETSC_SUCCESS);
1750: }

1752: /*@
1753:   TaoMonitorConstraintNorm - same as `TaoMonitorDefault()` except
1754:   it prints the norm of the constraint function.

1756:   Collective

1758:   Input Parameters:
1759: + tao - the `Tao` context
1760: - vf  - `PetscViewerAndFormat` context

1762:   Options Database Key:
1763: . -tao_monitor_constraint_norm - monitor the constraints

1765:   Level: advanced

1767: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1768: @*/
1769: PetscErrorCode TaoMonitorConstraintNorm(Tao tao, PetscViewerAndFormat *vf)
1770: {
1771:   PetscViewer viewer = vf->viewer;
1772:   PetscBool   isascii;
1773:   PetscInt    tabs;

1775:   PetscFunctionBegin;
1777:   if (vf->view_interval > 0 && tao->niter % vf->view_interval) PetscFunctionReturn(PETSC_SUCCESS);

1779:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1780:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1781:   if (isascii) {
1782:     PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1783:     PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1784:     PetscCall(PetscViewerASCIIPrintf(viewer, "iter = %" PetscInt_FMT ",", tao->niter));
1785:     PetscCall(PetscViewerASCIIPrintf(viewer, " Function value: %g,", (double)tao->fc));
1786:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g ", (double)tao->residual));
1787:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Constraint: %g \n", (double)tao->cnorm));
1788:     PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1789:   }
1790:   PetscCall(PetscViewerPopFormat(viewer));
1791:   PetscFunctionReturn(PETSC_SUCCESS);
1792: }

1794: /*@C
1795:   TaoMonitorSolution - Views the solution at each iteration of `TaoSolve()`

1797:   Collective

1799:   Input Parameters:
1800: + tao - the `Tao` context
1801: - vf  - `PetscViewerAndFormat` context

1803:   Options Database Key:
1804: . -tao_monitor_solution - view the solution

1806:   Level: advanced

1808: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1809: @*/
1810: PetscErrorCode TaoMonitorSolution(Tao tao, PetscViewerAndFormat *vf)
1811: {
1812:   PetscFunctionBegin;
1814:   if (vf->view_interval > 0 && tao->niter % vf->view_interval) PetscFunctionReturn(PETSC_SUCCESS);
1815:   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
1816:   PetscCall(VecView(tao->gradient, vf->viewer));
1817:   PetscCall(PetscViewerPopFormat(vf->viewer));
1818:   PetscFunctionReturn(PETSC_SUCCESS);
1819: }

1821: /*@C
1822:   TaoMonitorGradient - Views the gradient at each iteration of `TaoSolve()`

1824:   Collective

1826:   Input Parameters:
1827: + tao - the `Tao` context
1828: - vf  - `PetscViewerAndFormat` context

1830:   Options Database Key:
1831: . -tao_monitor_gradient - view the gradient at each iteration

1833:   Level: advanced

1835: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1836: @*/
1837: PetscErrorCode TaoMonitorGradient(Tao tao, PetscViewerAndFormat *vf)
1838: {
1839:   PetscFunctionBegin;
1841:   if (vf->view_interval > 0 && tao->niter % vf->view_interval) PetscFunctionReturn(PETSC_SUCCESS);
1842:   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
1843:   PetscCall(VecView(tao->gradient, vf->viewer));
1844:   PetscCall(PetscViewerPopFormat(vf->viewer));
1845:   PetscFunctionReturn(PETSC_SUCCESS);
1846: }

1848: /*@C
1849:   TaoMonitorStep - Views the step-direction at each iteration of `TaoSolve()`

1851:   Collective

1853:   Input Parameters:
1854: + tao - the `Tao` context
1855: - vf  - `PetscViewerAndFormat` context

1857:   Options Database Key:
1858: . -tao_monitor_step - view the step vector at each iteration

1860:   Level: advanced

1862: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1863: @*/
1864: PetscErrorCode TaoMonitorStep(Tao tao, PetscViewerAndFormat *vf)
1865: {
1866:   PetscFunctionBegin;
1868:   if (vf->view_interval > 0 && tao->niter % vf->view_interval) PetscFunctionReturn(PETSC_SUCCESS);
1869:   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
1870:   PetscCall(VecView(tao->stepdirection, vf->viewer));
1871:   PetscCall(PetscViewerPopFormat(vf->viewer));
1872:   PetscFunctionReturn(PETSC_SUCCESS);
1873: }

1875: /*@C
1876:   TaoMonitorSolutionDraw - Plots the solution at each iteration of `TaoSolve()`

1878:   Collective

1880:   Input Parameters:
1881: + tao - the `Tao` context
1882: - ctx - `TaoMonitorDraw` context

1884:   Options Database Key:
1885: . -tao_monitor_solution_draw - draw the solution at each iteration

1887:   Level: advanced

1889:   Note:
1890:   The context created by `TaoMonitorDrawCtxCreate()`, along with `TaoMonitorSolutionDraw()`, and `TaoMonitorDrawCtxDestroy()`
1891:   are passed to `TaoMonitorSet()` to monitor the solution graphically.

1893: .seealso: [](ch_tao), `Tao`, `TaoMonitorSolution()`, `TaoMonitorSet()`, `TaoMonitorGradientDraw()`, `TaoMonitorDrawCtxCreate()`,
1894:           `TaoMonitorDrawCtxDestroy()`
1895: @*/
1896: PetscErrorCode TaoMonitorSolutionDraw(Tao tao, PetscCtx ctx)
1897: {
1898:   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;

1900:   PetscFunctionBegin;
1902:   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(PETSC_SUCCESS);
1903:   PetscCall(VecView(tao->solution, ictx->viewer));
1904:   PetscFunctionReturn(PETSC_SUCCESS);
1905: }

1907: /*@C
1908:   TaoMonitorGradientDraw - Plots the gradient at each iteration of `TaoSolve()`

1910:   Collective

1912:   Input Parameters:
1913: + tao - the `Tao` context
1914: - ctx - `PetscViewer` context

1916:   Options Database Key:
1917: . -tao_monitor_gradient_draw - draw the gradient at each iteration

1919:   Level: advanced

1921: .seealso: [](ch_tao), `Tao`, `TaoMonitorGradient()`, `TaoMonitorSet()`, `TaoMonitorSolutionDraw()`
1922: @*/
1923: PetscErrorCode TaoMonitorGradientDraw(Tao tao, PetscCtx ctx)
1924: {
1925:   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;

1927:   PetscFunctionBegin;
1929:   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(PETSC_SUCCESS);
1930:   PetscCall(VecView(tao->gradient, ictx->viewer));
1931:   PetscFunctionReturn(PETSC_SUCCESS);
1932: }

1934: /*@C
1935:   TaoMonitorStepDraw - Plots the step direction at each iteration of `TaoSolve()`

1937:   Collective

1939:   Input Parameters:
1940: + tao - the `Tao` context
1941: - ctx - the `PetscViewer` context

1943:   Options Database Key:
1944: . -tao_monitor_step_draw - draw the step direction at each iteration

1946:   Level: advanced

1948: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorSolutionDraw`
1949: @*/
1950: PetscErrorCode TaoMonitorStepDraw(Tao tao, PetscCtx ctx)
1951: {
1952:   PetscViewer viewer = (PetscViewer)ctx;

1954:   PetscFunctionBegin;
1957:   PetscCall(VecView(tao->stepdirection, viewer));
1958:   PetscFunctionReturn(PETSC_SUCCESS);
1959: }

1961: /*@C
1962:   TaoMonitorResidual - Views the least-squares residual at each iteration of `TaoSolve()`

1964:   Collective

1966:   Input Parameters:
1967: + tao - the `Tao` context
1968: - vf  - `PetscViewerAndFormat` context

1970:   Options Database Key:
1971: . -tao_monitor_ls_residual - view the residual at each iteration

1973:   Level: advanced

1975: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1976: @*/
1977: PetscErrorCode TaoMonitorResidual(Tao tao, PetscViewerAndFormat *vf)
1978: {
1979:   PetscFunctionBegin;
1981:   if (vf->view_interval > 0 && tao->niter % vf->view_interval) PetscFunctionReturn(PETSC_SUCCESS);
1982:   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
1983:   PetscCall(VecView(tao->ls_res, vf->viewer));
1984:   PetscCall(PetscViewerPopFormat(vf->viewer));
1985:   PetscFunctionReturn(PETSC_SUCCESS);
1986: }

1988: /*@
1989:   TaoDefaultConvergenceTest - Determines whether the solver should continue iterating
1990:   or terminate.

1992:   Collective

1994:   Input Parameters:
1995: + tao   - the `Tao` context
1996: - dummy - unused dummy context

1998:   Level: developer

2000:   Notes:
2001:   This routine checks the residual in the optimality conditions, the
2002:   relative residual in the optimity conditions, the number of function
2003:   evaluations, and the function value to test convergence.  Some
2004:   solvers may use different convergence routines.

2006: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoGetConvergedReason()`, `TaoSetConvergedReason()`
2007: @*/
2008: PetscErrorCode TaoDefaultConvergenceTest(Tao tao, void *dummy)
2009: {
2010:   PetscInt           niter = tao->niter, nfuncs = PetscMax(tao->nfuncs, tao->nfuncgrads);
2011:   PetscInt           max_funcs = tao->max_funcs;
2012:   PetscReal          gnorm = tao->residual, gnorm0 = tao->gnorm0;
2013:   PetscReal          f = tao->fc, steptol = tao->steptol, trradius = tao->step;
2014:   PetscReal          gatol = tao->gatol, grtol = tao->grtol, gttol = tao->gttol;
2015:   PetscReal          catol = tao->catol, crtol = tao->crtol;
2016:   PetscReal          fmin = tao->fmin, cnorm = tao->cnorm;
2017:   TaoConvergedReason reason = tao->reason;

2019:   PetscFunctionBegin;
2021:   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);

2023:   if (PetscIsInfOrNanReal(f)) {
2024:     PetscCall(PetscInfo(tao, "Failed to converged, function value is infinity or NaN\n"));
2025:     reason = TAO_DIVERGED_NAN;
2026:   } else if (f <= fmin && cnorm <= catol) {
2027:     PetscCall(PetscInfo(tao, "Converged due to function value %g < minimum function value %g\n", (double)f, (double)fmin));
2028:     reason = TAO_CONVERGED_MINF;
2029:   } else if (gnorm <= gatol && cnorm <= catol) {
2030:     PetscCall(PetscInfo(tao, "Converged due to residual norm ||g(X)||=%g < %g\n", (double)gnorm, (double)gatol));
2031:     reason = TAO_CONVERGED_GATOL;
2032:   } else if (f != 0 && PetscAbsReal(gnorm / f) <= grtol && cnorm <= crtol) {
2033:     PetscCall(PetscInfo(tao, "Converged due to residual ||g(X)||/|f(X)| =%g < %g\n", (double)(gnorm / f), (double)grtol));
2034:     reason = TAO_CONVERGED_GRTOL;
2035:   } else if (gnorm0 != 0 && ((gttol == 0 && gnorm == 0) || gnorm / gnorm0 < gttol) && cnorm <= crtol) {
2036:     PetscCall(PetscInfo(tao, "Converged due to relative residual norm ||g(X)||/||g(X0)|| = %g < %g\n", (double)(gnorm / gnorm0), (double)gttol));
2037:     reason = TAO_CONVERGED_GTTOL;
2038:   } else if (max_funcs != PETSC_UNLIMITED && nfuncs > max_funcs) {
2039:     PetscCall(PetscInfo(tao, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", nfuncs, max_funcs));
2040:     reason = TAO_DIVERGED_MAXFCN;
2041:   } else if (tao->lsflag != 0) {
2042:     PetscCall(PetscInfo(tao, "Tao Line Search failure.\n"));
2043:     reason = TAO_DIVERGED_LS_FAILURE;
2044:   } else if (trradius < steptol && niter > 0) {
2045:     PetscCall(PetscInfo(tao, "Trust region/step size too small: %g < %g\n", (double)trradius, (double)steptol));
2046:     reason = TAO_CONVERGED_STEPTOL;
2047:   } else if (niter >= tao->max_it) {
2048:     PetscCall(PetscInfo(tao, "Exceeded maximum number of iterations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", niter, tao->max_it));
2049:     reason = TAO_DIVERGED_MAXITS;
2050:   } else {
2051:     reason = TAO_CONTINUE_ITERATING;
2052:   }
2053:   tao->reason = reason;
2054:   PetscFunctionReturn(PETSC_SUCCESS);
2055: }

2057: /*@
2058:   TaoSetOptionsPrefix - Sets the prefix used for searching for all
2059:   Tao options in the database.

2061:   Logically Collective

2063:   Input Parameters:
2064: + tao - the `Tao` context
2065: - p   - the prefix string to prepend to all Tao option requests

2067:   Level: advanced

2069:   Notes:
2070:   A hyphen (-) must NOT be given at the beginning of the prefix name.
2071:   The first character of all runtime options is AUTOMATICALLY the hyphen.

2073:   For example, to distinguish between the runtime options for two
2074:   different Tao solvers, one could call
2075: .vb
2076:       TaoSetOptionsPrefix(tao1,"sys1_")
2077:       TaoSetOptionsPrefix(tao2,"sys2_")
2078: .ve

2080:   This would enable use of different options for each system, such as
2081: .vb
2082:       -sys1_tao_method blmvm -sys1_tao_grtol 1.e-3
2083:       -sys2_tao_method lmvm  -sys2_tao_grtol 1.e-4
2084: .ve

2086: .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoAppendOptionsPrefix()`, `TaoGetOptionsPrefix()`
2087: @*/
2088: PetscErrorCode TaoSetOptionsPrefix(Tao tao, const char p[])
2089: {
2090:   PetscFunctionBegin;
2092:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tao, p));
2093:   if (tao->linesearch) PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, p));
2094:   if (tao->ksp) PetscCall(KSPSetOptionsPrefix(tao->ksp, p));
2095:   PetscFunctionReturn(PETSC_SUCCESS);
2096: }

2098: /*@
2099:   TaoAppendOptionsPrefix - Appends to the prefix used for searching for all Tao options in the database.

2101:   Logically Collective

2103:   Input Parameters:
2104: + tao - the `Tao` solver context
2105: - p   - the prefix string to prepend to all `Tao` option requests

2107:   Level: advanced

2109:   Note:
2110:   A hyphen (-) must NOT be given at the beginning of the prefix name.
2111:   The first character of all runtime options is automatically the hyphen.

2113: .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoSetOptionsPrefix()`, `TaoGetOptionsPrefix()`
2114: @*/
2115: PetscErrorCode TaoAppendOptionsPrefix(Tao tao, const char p[])
2116: {
2117:   PetscFunctionBegin;
2119:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao, p));
2120:   if (tao->linesearch) PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao->linesearch, p));
2121:   if (tao->ksp) PetscCall(KSPAppendOptionsPrefix(tao->ksp, p));
2122:   PetscFunctionReturn(PETSC_SUCCESS);
2123: }

2125: /*@
2126:   TaoGetOptionsPrefix - Gets the prefix used for searching for all
2127:   Tao options in the database

2129:   Not Collective

2131:   Input Parameter:
2132: . tao - the `Tao` context

2134:   Output Parameter:
2135: . p - pointer to the prefix string used is returned

2137:   Level: advanced

2139: .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoSetOptionsPrefix()`, `TaoAppendOptionsPrefix()`
2140: @*/
2141: PetscErrorCode TaoGetOptionsPrefix(Tao tao, const char *p[])
2142: {
2143:   PetscFunctionBegin;
2145:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)tao, p));
2146:   PetscFunctionReturn(PETSC_SUCCESS);
2147: }

2149: /*@
2150:   TaoSetType - Sets the `TaoType` for the minimization solver.

2152:   Collective

2154:   Input Parameters:
2155: + tao  - the `Tao` solver context
2156: - type - a known method

2158:   Options Database Key:
2159: . -tao_type <type> - Sets the method; use -help for a list
2160:    of available methods (for instance, "-tao_type lmvm" or "-tao_type tron")

2162:   Level: intermediate

2164: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoGetType()`, `TaoType`
2165: @*/
2166: PetscErrorCode TaoSetType(Tao tao, TaoType type)
2167: {
2168:   PetscErrorCode (*create_xxx)(Tao);
2169:   PetscBool issame;

2171:   PetscFunctionBegin;

2174:   PetscCall(PetscObjectTypeCompare((PetscObject)tao, type, &issame));
2175:   if (issame) PetscFunctionReturn(PETSC_SUCCESS);

2177:   PetscCall(PetscFunctionListFind(TaoList, type, &create_xxx));
2178:   PetscCheck(create_xxx, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Tao type %s", type);

2180:   /* Destroy the existing solver information */
2181:   PetscTryTypeMethod(tao, destroy);
2182:   PetscCall(KSPDestroy(&tao->ksp));
2183:   PetscCall(TaoLineSearchDestroy(&tao->linesearch));

2185:   /* Reinitialize type-specific function pointers in TaoOps structure */
2186:   tao->ops->setup          = NULL;
2187:   tao->ops->computedual    = NULL;
2188:   tao->ops->solve          = NULL;
2189:   tao->ops->view           = NULL;
2190:   tao->ops->setfromoptions = NULL;
2191:   tao->ops->destroy        = NULL;

2193:   tao->setupcalled = PETSC_FALSE;

2195:   PetscCall((*create_xxx)(tao));
2196:   PetscCall(PetscObjectChangeTypeName((PetscObject)tao, type));
2197:   PetscFunctionReturn(PETSC_SUCCESS);
2198: }

2200: /*@C
2201:   TaoRegister - Adds a method to the Tao package for minimization.

2203:   Not Collective, No Fortran Support

2205:   Input Parameters:
2206: + sname - name of a new user-defined solver
2207: - func  - routine to create `TaoType` specific method context

2209:   Calling sequence of `func`:
2210: . tao - the `Tao` object to be created

2212:   Example Usage:
2213: .vb
2214:    TaoRegister("my_solver", MySolverCreate);
2215: .ve

2217:   Then, your solver can be chosen with the procedural interface via
2218: .vb
2219:   TaoSetType(tao, "my_solver")
2220: .ve
2221:   or at runtime via the option
2222: .vb
2223:   -tao_type my_solver
2224: .ve

2226:   Level: advanced

2228:   Note:
2229:   `TaoRegister()` may be called multiple times to add several user-defined solvers.

2231: .seealso: [](ch_tao), `Tao`, `TaoSetType()`, `TaoRegisterAll()`, `TaoRegisterDestroy()`
2232: @*/
2233: PetscErrorCode TaoRegister(const char sname[], PetscErrorCode (*func)(Tao tao))
2234: {
2235:   PetscFunctionBegin;
2236:   PetscCall(TaoInitializePackage());
2237:   PetscCall(PetscFunctionListAdd(&TaoList, sname, func));
2238:   PetscFunctionReturn(PETSC_SUCCESS);
2239: }

2241: /*@C
2242:   TaoRegisterDestroy - Frees the list of minimization solvers that were
2243:   registered by `TaoRegister()`.

2245:   Not Collective

2247:   Level: advanced

2249: .seealso: [](ch_tao), `Tao`, `TaoRegisterAll()`, `TaoRegister()`
2250: @*/
2251: PetscErrorCode TaoRegisterDestroy(void)
2252: {
2253:   PetscFunctionBegin;
2254:   PetscCall(PetscFunctionListDestroy(&TaoList));
2255:   TaoRegisterAllCalled = PETSC_FALSE;
2256:   PetscFunctionReturn(PETSC_SUCCESS);
2257: }

2259: /*@
2260:   TaoGetIterationNumber - Gets the number of `TaoSolve()` iterations completed
2261:   at this time.

2263:   Not Collective

2265:   Input Parameter:
2266: . tao - the `Tao` context

2268:   Output Parameter:
2269: . iter - iteration number

2271:   Notes:
2272:   For example, during the computation of iteration 2 this would return 1.

2274:   Level: intermediate

2276: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`, `TaoGetResidualNorm()`, `TaoGetObjective()`
2277: @*/
2278: PetscErrorCode TaoGetIterationNumber(Tao tao, PetscInt *iter)
2279: {
2280:   PetscFunctionBegin;
2282:   PetscAssertPointer(iter, 2);
2283:   *iter = tao->niter;
2284:   PetscFunctionReturn(PETSC_SUCCESS);
2285: }

2287: /*@
2288:   TaoGetResidualNorm - Gets the current value of the norm of the residual (gradient)
2289:   at this time.

2291:   Not Collective

2293:   Input Parameter:
2294: . tao - the `Tao` context

2296:   Output Parameter:
2297: . value - the current value

2299:   Level: intermediate

2301:   Developer Notes:
2302:   This is the 2-norm of the residual, we cannot use `TaoGetGradientNorm()` because that has
2303:   a different meaning. For some reason `Tao` sometimes calls the gradient the residual.

2305: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`, `TaoGetIterationNumber()`, `TaoGetObjective()`
2306: @*/
2307: PetscErrorCode TaoGetResidualNorm(Tao tao, PetscReal *value)
2308: {
2309:   PetscFunctionBegin;
2311:   PetscAssertPointer(value, 2);
2312:   *value = tao->residual;
2313:   PetscFunctionReturn(PETSC_SUCCESS);
2314: }

2316: /*@
2317:   TaoSetIterationNumber - Sets the current iteration number.

2319:   Logically Collective

2321:   Input Parameters:
2322: + tao  - the `Tao` context
2323: - iter - iteration number

2325:   Level: developer

2327: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
2328: @*/
2329: PetscErrorCode TaoSetIterationNumber(Tao tao, PetscInt iter)
2330: {
2331:   PetscFunctionBegin;
2334:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2335:   tao->niter = iter;
2336:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2337:   PetscFunctionReturn(PETSC_SUCCESS);
2338: }

2340: /*@
2341:   TaoGetTotalIterationNumber - Gets the total number of `TaoSolve()` iterations
2342:   completed. This number keeps accumulating if multiple solves
2343:   are called with the `Tao` object.

2345:   Not Collective

2347:   Input Parameter:
2348: . tao - the `Tao` context

2350:   Output Parameter:
2351: . iter - number of iterations

2353:   Level: intermediate

2355:   Note:
2356:   The total iteration count is updated after each solve, if there is a current
2357:   `TaoSolve()` in progress then those iterations are not included in the count

2359: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
2360: @*/
2361: PetscErrorCode TaoGetTotalIterationNumber(Tao tao, PetscInt *iter)
2362: {
2363:   PetscFunctionBegin;
2365:   PetscAssertPointer(iter, 2);
2366:   *iter = tao->ntotalits;
2367:   PetscFunctionReturn(PETSC_SUCCESS);
2368: }

2370: /*@
2371:   TaoSetTotalIterationNumber - Sets the current total iteration number.

2373:   Logically Collective

2375:   Input Parameters:
2376: + tao  - the `Tao` context
2377: - iter - the iteration number

2379:   Level: developer

2381: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
2382: @*/
2383: PetscErrorCode TaoSetTotalIterationNumber(Tao tao, PetscInt iter)
2384: {
2385:   PetscFunctionBegin;
2388:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2389:   tao->ntotalits = iter;
2390:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2391:   PetscFunctionReturn(PETSC_SUCCESS);
2392: }

2394: /*@
2395:   TaoSetConvergedReason - Sets the termination flag on a `Tao` object

2397:   Logically Collective

2399:   Input Parameters:
2400: + tao    - the `Tao` context
2401: - reason - the `TaoConvergedReason`

2403:   Level: intermediate

2405: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`
2406: @*/
2407: PetscErrorCode TaoSetConvergedReason(Tao tao, TaoConvergedReason reason)
2408: {
2409:   PetscFunctionBegin;
2412:   tao->reason = reason;
2413:   PetscFunctionReturn(PETSC_SUCCESS);
2414: }

2416: /*@
2417:   TaoGetConvergedReason - Gets the reason the `TaoSolve()` was stopped.

2419:   Not Collective

2421:   Input Parameter:
2422: . tao - the `Tao` solver context

2424:   Output Parameter:
2425: . reason - value of `TaoConvergedReason`

2427:   Level: intermediate

2429: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetConvergenceTest()`, `TaoSetTolerances()`
2430: @*/
2431: PetscErrorCode TaoGetConvergedReason(Tao tao, TaoConvergedReason *reason)
2432: {
2433:   PetscFunctionBegin;
2435:   PetscAssertPointer(reason, 2);
2436:   *reason = tao->reason;
2437:   PetscFunctionReturn(PETSC_SUCCESS);
2438: }

2440: /*@
2441:   TaoGetSolutionStatus - Get the current iterate, objective value,
2442:   residual, infeasibility, and termination from a `Tao` object

2444:   Not Collective

2446:   Input Parameter:
2447: . tao - the `Tao` context

2449:   Output Parameters:
2450: + its    - the current iterate number (>=0)
2451: . f      - the current function value
2452: . gnorm  - the square of the gradient norm, duality gap, or other measure indicating distance from optimality.
2453: . cnorm  - the infeasibility of the current solution with regard to the constraints.
2454: . xdiff  - the step length or trust region radius of the most recent iterate.
2455: - reason - The termination reason, which can equal `TAO_CONTINUE_ITERATING`

2457:   Level: intermediate

2459:   Notes:
2460:   Tao returns the values set by the solvers in the routine `TaoMonitor()`.

2462:   If any of the output arguments are set to `NULL`, no corresponding value will be returned.

2464: .seealso: [](ch_tao), `TaoMonitor()`, `TaoGetConvergedReason()`
2465: @*/
2466: PetscErrorCode TaoGetSolutionStatus(Tao tao, PetscInt *its, PetscReal *f, PetscReal *gnorm, PetscReal *cnorm, PetscReal *xdiff, TaoConvergedReason *reason)
2467: {
2468:   PetscFunctionBegin;
2470:   if (its) *its = tao->niter;
2471:   if (f) *f = tao->fc;
2472:   if (gnorm) *gnorm = tao->residual;
2473:   if (cnorm) *cnorm = tao->cnorm;
2474:   if (reason) *reason = tao->reason;
2475:   if (xdiff) *xdiff = tao->step;
2476:   PetscFunctionReturn(PETSC_SUCCESS);
2477: }

2479: /*@
2480:   TaoGetType - Gets the current `TaoType` being used in the `Tao` object

2482:   Not Collective

2484:   Input Parameter:
2485: . tao - the `Tao` solver context

2487:   Output Parameter:
2488: . type - the `TaoType`

2490:   Level: intermediate

2492: .seealso: [](ch_tao), `Tao`, `TaoType`, `TaoSetType()`
2493: @*/
2494: PetscErrorCode TaoGetType(Tao tao, TaoType *type)
2495: {
2496:   PetscFunctionBegin;
2498:   PetscAssertPointer(type, 2);
2499:   *type = ((PetscObject)tao)->type_name;
2500:   PetscFunctionReturn(PETSC_SUCCESS);
2501: }

2503: /*@C
2504:   TaoMonitor - Monitor the solver and the current solution.  This
2505:   routine will record the iteration number and residual statistics,
2506:   and call any monitors specified by the user.

2508:   Input Parameters:
2509: + tao        - the `Tao` context
2510: . its        - the current iterate number (>=0)
2511: . f          - the current objective function value
2512: . res        - the gradient norm, square root of the duality gap, or other measure indicating distance from optimality.  This measure will be recorded and
2513:           used for some termination tests.
2514: . cnorm      - the infeasibility of the current solution with regard to the constraints.
2515: - steplength - multiple of the step direction added to the previous iterate.

2517:   Options Database Key:
2518: . -tao_monitor - Use the default monitor, which prints statistics to standard output

2520:   Level: developer

2522: .seealso: [](ch_tao), `Tao`, `TaoGetConvergedReason()`, `TaoMonitorDefault()`, `TaoMonitorSet()`
2523: @*/
2524: PetscErrorCode TaoMonitor(Tao tao, PetscInt its, PetscReal f, PetscReal res, PetscReal cnorm, PetscReal steplength)
2525: {
2526:   PetscInt i;

2528:   PetscFunctionBegin;
2530:   tao->fc       = f;
2531:   tao->residual = res;
2532:   tao->cnorm    = cnorm;
2533:   tao->step     = steplength;
2534:   if (!its) {
2535:     tao->cnorm0 = cnorm;
2536:     tao->gnorm0 = res;
2537:   }
2538:   PetscCall(VecLockReadPush(tao->solution));
2539:   for (i = 0; i < tao->numbermonitors; i++) PetscCall((*tao->monitor[i])(tao, tao->monitorcontext[i]));
2540:   PetscCall(VecLockReadPop(tao->solution));
2541:   PetscFunctionReturn(PETSC_SUCCESS);
2542: }

2544: /*@
2545:   TaoSetConvergenceHistory - Sets the array used to hold the convergence history.

2547:   Logically Collective

2549:   Input Parameters:
2550: + tao   - the `Tao` solver context
2551: . obj   - array to hold objective value history
2552: . resid - array to hold residual history
2553: . cnorm - array to hold constraint violation history
2554: . lits  - integer array holds the number of linear iterations for each Tao iteration
2555: . na    - size of `obj`, `resid`, and `cnorm`
2556: - reset - `PETSC_TRUE` indicates each new minimization resets the history counter to zero,
2557:            else it continues storing new values for new minimizations after the old ones

2559:   Level: intermediate

2561:   Notes:
2562:   If set, `Tao` will fill the given arrays with the indicated
2563:   information at each iteration.  If 'obj','resid','cnorm','lits' are
2564:   *all* `NULL` then space (using size `na`, or 1000 if `na` is `PETSC_DECIDE`) is allocated for the history.
2565:   If not all are `NULL`, then only the non-`NULL` information categories
2566:   will be stored, the others will be ignored.

2568:   Any convergence information after iteration number 'na' will not be stored.

2570:   This routine is useful, e.g., when running a code for purposes
2571:   of accurate performance monitoring, when no I/O should be done
2572:   during the section of code that is being timed.

2574: .seealso: [](ch_tao), `TaoGetConvergenceHistory()`
2575: @*/
2576: PetscErrorCode TaoSetConvergenceHistory(Tao tao, PetscReal obj[], PetscReal resid[], PetscReal cnorm[], PetscInt lits[], PetscInt na, PetscBool reset)
2577: {
2578:   PetscFunctionBegin;
2580:   if (obj) PetscAssertPointer(obj, 2);
2581:   if (resid) PetscAssertPointer(resid, 3);
2582:   if (cnorm) PetscAssertPointer(cnorm, 4);
2583:   if (lits) PetscAssertPointer(lits, 5);

2585:   if (na == PETSC_DECIDE || na == PETSC_CURRENT) na = 1000;
2586:   if (!obj && !resid && !cnorm && !lits) {
2587:     PetscCall(PetscCalloc4(na, &obj, na, &resid, na, &cnorm, na, &lits));
2588:     tao->hist_malloc = PETSC_TRUE;
2589:   }

2591:   tao->hist_obj   = obj;
2592:   tao->hist_resid = resid;
2593:   tao->hist_cnorm = cnorm;
2594:   tao->hist_lits  = lits;
2595:   tao->hist_max   = na;
2596:   tao->hist_reset = reset;
2597:   tao->hist_len   = 0;
2598:   PetscFunctionReturn(PETSC_SUCCESS);
2599: }

2601: /*@C
2602:   TaoGetConvergenceHistory - Gets the arrays used that hold the convergence history.

2604:   Collective

2606:   Input Parameter:
2607: . tao - the `Tao` context

2609:   Output Parameters:
2610: + obj   - array used to hold objective value history
2611: . resid - array used to hold residual history
2612: . cnorm - array used to hold constraint violation history
2613: . lits  - integer array used to hold linear solver iteration count
2614: - nhist - size of `obj`, `resid`, `cnorm`, and `lits`

2616:   Level: advanced

2618:   Notes:
2619:   This routine must be preceded by calls to `TaoSetConvergenceHistory()`
2620:   and `TaoSolve()`, otherwise it returns useless information.

2622:   This routine is useful, e.g., when running a code for purposes
2623:   of accurate performance monitoring, when no I/O should be done
2624:   during the section of code that is being timed.

2626:   Fortran Notes:
2627:   The calling sequence is
2628: .vb
2629:    call TaoGetConvergenceHistory(Tao tao, PetscInt nhist, PetscErrorCode ierr)
2630: .ve
2631:   In other words this gets the current number of entries in the history. Access the history through the array you passed to `TaoSetConvergenceHistory()`

2633: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetConvergenceHistory()`
2634: @*/
2635: PetscErrorCode TaoGetConvergenceHistory(Tao tao, PetscReal **obj, PetscReal **resid, PetscReal **cnorm, PetscInt **lits, PetscInt *nhist)
2636: {
2637:   PetscFunctionBegin;
2639:   if (obj) *obj = tao->hist_obj;
2640:   if (cnorm) *cnorm = tao->hist_cnorm;
2641:   if (resid) *resid = tao->hist_resid;
2642:   if (lits) *lits = tao->hist_lits;
2643:   if (nhist) *nhist = tao->hist_len;
2644:   PetscFunctionReturn(PETSC_SUCCESS);
2645: }

2647: /*@
2648:   TaoSetApplicationContext - Sets the optional user-defined context for a `Tao` solver that can be accessed later, for example in the
2649:   `Tao` callback functions with `TaoGetApplicationContext()`

2651:   Logically Collective

2653:   Input Parameters:
2654: + tao - the `Tao` context
2655: - ctx - the user context

2657:   Level: intermediate

2659:   Fortran Note:
2660:   This only works when `ctx` is a Fortran derived type (it cannot be a `PetscObject`), we recommend writing a Fortran interface definition for this
2661:   function that tells the Fortran compiler the derived data type that is passed in as the `ctx` argument. See `TaoGetApplicationContext()` for
2662:   an example.

2664: .seealso: [](ch_tao), `Tao`, `TaoGetApplicationContext()`
2665: @*/
2666: PetscErrorCode TaoSetApplicationContext(Tao tao, PetscCtx ctx)
2667: {
2668:   PetscFunctionBegin;
2670:   tao->ctx = ctx;
2671:   PetscFunctionReturn(PETSC_SUCCESS);
2672: }

2674: /*@
2675:   TaoGetApplicationContext - Gets the user-defined context for a `Tao` solver provided with `TaoSetApplicationContext()`

2677:   Not Collective

2679:   Input Parameter:
2680: . tao - the `Tao` context

2682:   Output Parameter:
2683: . ctx - a pointer to the user context

2685:   Level: intermediate

2687:   Fortran Note:
2688:   This only works when the context is a Fortran derived type or a `PetscObject`. Define `ctx` with
2689: .vb
2690:   type(tUsertype), pointer :: ctx
2691: .ve

2693: .seealso: [](ch_tao), `Tao`, `TaoSetApplicationContext()`
2694: @*/
2695: PetscErrorCode TaoGetApplicationContext(Tao tao, PetscCtxRt ctx)
2696: {
2697:   PetscFunctionBegin;
2699:   PetscAssertPointer(ctx, 2);
2700:   *(void **)ctx = tao->ctx;
2701:   PetscFunctionReturn(PETSC_SUCCESS);
2702: }

2704: /*@
2705:   TaoSetGradientNorm - Sets the matrix used to define the norm that measures the size of the gradient in some of the `Tao` algorithms

2707:   Collective

2709:   Input Parameters:
2710: + tao - the `Tao` context
2711: - M   - matrix that defines the norm

2713:   Level: beginner

2715: .seealso: [](ch_tao), `Tao`, `TaoGetGradientNorm()`, `TaoGradientNorm()`
2716: @*/
2717: PetscErrorCode TaoSetGradientNorm(Tao tao, Mat M)
2718: {
2719:   PetscFunctionBegin;
2722:   PetscCall(PetscObjectReference((PetscObject)M));
2723:   PetscCall(MatDestroy(&tao->gradient_norm));
2724:   PetscCall(VecDestroy(&tao->gradient_norm_tmp));
2725:   tao->gradient_norm = M;
2726:   PetscCall(MatCreateVecs(M, NULL, &tao->gradient_norm_tmp));
2727:   PetscFunctionReturn(PETSC_SUCCESS);
2728: }

2730: /*@
2731:   TaoGetGradientNorm - Returns the matrix used to define the norm used for measuring the size of the gradient in some of the `Tao` algorithms

2733:   Not Collective

2735:   Input Parameter:
2736: . tao - the `Tao` context

2738:   Output Parameter:
2739: . M - gradient norm

2741:   Level: beginner

2743: .seealso: [](ch_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGradientNorm()`
2744: @*/
2745: PetscErrorCode TaoGetGradientNorm(Tao tao, Mat *M)
2746: {
2747:   PetscFunctionBegin;
2749:   PetscAssertPointer(M, 2);
2750:   *M = tao->gradient_norm;
2751:   PetscFunctionReturn(PETSC_SUCCESS);
2752: }

2754: /*@
2755:   TaoGradientNorm - Compute the norm using the `NormType`, the user has selected

2757:   Collective

2759:   Input Parameters:
2760: + tao      - the `Tao` context
2761: . gradient - the gradient
2762: - type     - the norm type

2764:   Output Parameter:
2765: . gnorm - the gradient norm

2767:   Level: advanced

2769:   Note:
2770:   If `TaoSetGradientNorm()` has been set and `type` is `NORM_2` then the norm provided with `TaoSetGradientNorm()` is used.

2772:   Developer Notes:
2773:   Should be named `TaoComputeGradientNorm()`.

2775:   The usage is a bit confusing, with `TaoSetGradientNorm()` plus `NORM_2` resulting in the computation of the user provided
2776:   norm, perhaps a refactorization is in order.

2778: .seealso: [](ch_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGetGradientNorm()`
2779: @*/
2780: PetscErrorCode TaoGradientNorm(Tao tao, Vec gradient, NormType type, PetscReal *gnorm)
2781: {
2782:   PetscFunctionBegin;
2786:   PetscAssertPointer(gnorm, 4);
2787:   if (tao->gradient_norm) {
2788:     PetscScalar gnorms;

2790:     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.");
2791:     PetscCall(MatMult(tao->gradient_norm, gradient, tao->gradient_norm_tmp));
2792:     PetscCall(VecDot(gradient, tao->gradient_norm_tmp, &gnorms));
2793:     *gnorm = PetscRealPart(PetscSqrtScalar(gnorms));
2794:   } else {
2795:     PetscCall(VecNorm(gradient, type, gnorm));
2796:   }
2797:   PetscFunctionReturn(PETSC_SUCCESS);
2798: }

2800: /*@C
2801:   TaoMonitorDrawCtxCreate - Creates the monitor context for `TaoMonitorSolutionDraw()`

2803:   Collective

2805:   Input Parameters:
2806: + comm     - the communicator to share the context
2807: . host     - the name of the X Windows host that will display the monitor
2808: . label    - the label to put at the top of the display window
2809: . x        - the horizontal coordinate of the lower left corner of the window to open
2810: . y        - the vertical coordinate of the lower left corner of the window to open
2811: . m        - the width of the window
2812: . n        - the height of the window
2813: - howoften - how many `Tao` iterations between displaying the monitor information

2815:   Output Parameter:
2816: . ctx - the monitor context

2818:   Options Database Keys:
2819: + -tao_monitor_solution_draw - use `TaoMonitorSolutionDraw()` to monitor the solution
2820: - -tao_draw_solution_initial - show initial guess as well as current solution

2822:   Level: intermediate

2824:   Note:
2825:   The context this creates, along with `TaoMonitorSolutionDraw()`, and `TaoMonitorDrawCtxDestroy()`
2826:   are passed to `TaoMonitorSet()`.

2828: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorDrawCtx()`
2829: @*/
2830: PetscErrorCode TaoMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TaoMonitorDrawCtx *ctx)
2831: {
2832:   PetscFunctionBegin;
2833:   PetscCall(PetscNew(ctx));
2834:   PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
2835:   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
2836:   (*ctx)->howoften = howoften;
2837:   PetscFunctionReturn(PETSC_SUCCESS);
2838: }

2840: /*@C
2841:   TaoMonitorDrawCtxDestroy - Destroys the monitor context for `TaoMonitorSolutionDraw()`

2843:   Collective

2845:   Input Parameter:
2846: . ictx - the monitor context

2848:   Level: intermediate

2850:   Note:
2851:   This is passed to `TaoMonitorSet()` as the final argument, along with `TaoMonitorSolutionDraw()`, and the context
2852:   obtained with `TaoMonitorDrawCtxCreate()`.

2854: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorSolutionDraw()`
2855: @*/
2856: PetscErrorCode TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx *ictx)
2857: {
2858:   PetscFunctionBegin;
2859:   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
2860:   PetscCall(PetscFree(*ictx));
2861:   PetscFunctionReturn(PETSC_SUCCESS);
2862: }