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, void *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, void *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);
187:     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)tao)->tablevel));
188:     if (tao->reason > 0) {
189:       PetscCall(PetscViewerASCIIPrintf(viewer, "  TAO %s solve converged due to %s iterations %" PetscInt_FMT "\n", ((PetscObject)tao)->prefix ? ((PetscObject)tao)->prefix : "", TaoConvergedReasons[tao->reason], tao->niter));
190:     } else {
191:       PetscCall(PetscViewerASCIIPrintf(viewer, "  TAO %s solve did not converge due to %s iteration %" PetscInt_FMT "\n", ((PetscObject)tao)->prefix ? ((PetscObject)tao)->prefix : "", TaoConvergedReasons[tao->reason], tao->niter));
192:     }
193:     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)tao)->tablevel));
194:   }
195:   PetscCall(TaoViewFromOptions(tao, NULL, "-tao_view"));
196:   PetscFunctionReturn(PETSC_SUCCESS);
197: }

199: /*@
200:   TaoSetUp - Sets up the internal data structures for the later use
201:   of a Tao solver

203:   Collective

205:   Input Parameter:
206: . tao - the `Tao` context

208:   Level: advanced

210:   Note:
211:   The user will not need to explicitly call `TaoSetUp()`, as it will
212:   automatically be called in `TaoSolve()`.  However, if the user
213:   desires to call it explicitly, it should come after `TaoCreate()`
214:   and any TaoSetSomething() routines, but before `TaoSolve()`.

216: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
217: @*/
218: PetscErrorCode TaoSetUp(Tao tao)
219: {
220:   PetscFunctionBegin;
222:   if (tao->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
223:   PetscCall(TaoSetUpEW_Private(tao));
224:   PetscCheck(tao->solution, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoSetSolution");
225:   PetscTryTypeMethod(tao, setup);
226:   tao->setupcalled = PETSC_TRUE;
227:   PetscFunctionReturn(PETSC_SUCCESS);
228: }

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

233:   Collective

235:   Input Parameter:
236: . tao - the `Tao` context

238:   Level: beginner

240: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
241: @*/
242: PetscErrorCode TaoDestroy(Tao *tao)
243: {
244:   PetscFunctionBegin;
245:   if (!*tao) PetscFunctionReturn(PETSC_SUCCESS);
247:   if (--((PetscObject)*tao)->refct > 0) {
248:     *tao = NULL;
249:     PetscFunctionReturn(PETSC_SUCCESS);
250:   }

252:   PetscTryTypeMethod(*tao, destroy);
253:   PetscCall(KSPDestroy(&(*tao)->ksp));
254:   PetscCall(SNESDestroy(&(*tao)->snes_ewdummy));
255:   PetscCall(TaoLineSearchDestroy(&(*tao)->linesearch));

257:   if ((*tao)->ops->convergencedestroy) {
258:     PetscCall((*(*tao)->ops->convergencedestroy)((*tao)->cnvP));
259:     if ((*tao)->jacobian_state_inv) PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
260:   }
261:   PetscCall(VecDestroy(&(*tao)->solution));
262:   PetscCall(VecDestroy(&(*tao)->gradient));
263:   PetscCall(VecDestroy(&(*tao)->ls_res));

265:   if ((*tao)->gradient_norm) {
266:     PetscCall(PetscObjectDereference((PetscObject)(*tao)->gradient_norm));
267:     PetscCall(VecDestroy(&(*tao)->gradient_norm_tmp));
268:   }

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

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

311:   Logically Collective

313:   Input Parameters:
314: + tao  - Tao context
315: - flag - `PETSC_TRUE` or `PETSC_FALSE`

317:   Level: advanced

319:   Note:
320:   See `SNESKSPSetUseEW()` for customization details.

322: .seealso: [](ch_tao), `Tao`, `SNESKSPSetUseEW()`
323: @*/
324: PetscErrorCode TaoKSPSetUseEW(Tao tao, PetscBool flag)
325: {
326:   PetscFunctionBegin;
329:   tao->ksp_ewconv = flag;
330:   PetscFunctionReturn(PETSC_SUCCESS);
331: }

333: /*@
334:   TaoSetFromOptions - Sets various Tao parameters from the options database

336:   Collective

338:   Input Parameter:
339: . tao - the `Tao` solver context

341:   Options Database Keys:
342: + -tao_type <type>             - The algorithm that Tao uses (lmvm, nls, etc.)
343: . -tao_gatol <gatol>           - absolute error tolerance for ||gradient||
344: . -tao_grtol <grtol>           - relative error tolerance for ||gradient||
345: . -tao_gttol <gttol>           - reduction of ||gradient|| relative to initial gradient
346: . -tao_max_it <max>            - sets maximum number of iterations
347: . -tao_max_funcs <max>         - sets maximum number of function evaluations
348: . -tao_fmin <fmin>             - stop if function value reaches fmin
349: . -tao_steptol <tol>           - stop if trust region radius less than <tol>
350: . -tao_trust0 <t>              - initial trust region radius
351: . -tao_view_solution           - view the solution at the end of the optimization process
352: . -tao_monitor                 - prints function value and residual norm at each iteration
353: . -tao_monitor_short           - same as `-tao_monitor`, but truncates very small values
354: . -tao_monitor_constraint_norm - prints objective value, gradient, and constraint norm at each iteration
355: . -tao_monitor_globalization   - prints information about the globalization at each iteration
356: . -tao_monitor_solution        - prints solution vector at each iteration
357: . -tao_monitor_ls_residual     - prints least-squares residual vector at each iteration
358: . -tao_monitor_step            - prints step vector at each iteration
359: . -tao_monitor_gradient        - prints gradient vector at each iteration
360: . -tao_monitor_solution_draw   - graphically view solution vector at each iteration
361: . -tao_monitor_step_draw       - graphically view step vector at each iteration
362: . -tao_monitor_gradient_draw   - graphically view gradient at each iteration
363: . -tao_monitor_cancel          - cancels all monitors (except those set with command line)
364: . -tao_fd_gradient             - use gradient computed with finite differences
365: . -tao_fd_hessian              - use hessian computed with finite differences
366: . -tao_mf_hessian              - use matrix-free Hessian computed with finite differences
367: . -tao_view                    - prints information about the Tao after solving
368: - -tao_converged_reason        - prints the reason Tao stopped iterating

370:   Level: beginner

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

376: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
377: @*/
378: PetscErrorCode TaoSetFromOptions(Tao tao)
379: {
380:   TaoType     default_type = TAOLMVM;
381:   char        type[256], monfilename[PETSC_MAX_PATH_LEN];
382:   PetscViewer monviewer;
383:   PetscBool   flg, found;
384:   MPI_Comm    comm;
385:   PetscReal   catol, crtol, gatol, grtol, gttol;

387:   PetscFunctionBegin;
389:   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));

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

393:   PetscObjectOptionsBegin((PetscObject)tao);
394:   /* Check for type from options */
395:   PetscCall(PetscOptionsFList("-tao_type", "Tao Solver type", "TaoSetType", TaoList, default_type, type, 256, &flg));
396:   if (flg) {
397:     PetscCall(TaoSetType(tao, type));
398:   } else if (!((PetscObject)tao)->type_name) {
399:     PetscCall(TaoSetType(tao, default_type));
400:   }

402:   /* Tao solvers do not set the prefix, set it here if not yet done
403:      We do it after SetType since solver may have been changed */
404:   if (tao->linesearch) {
405:     const char *prefix;
406:     PetscCall(TaoLineSearchGetOptionsPrefix(tao->linesearch, &prefix));
407:     if (!prefix) PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, ((PetscObject)tao)->prefix));
408:   }

410:   catol = tao->catol;
411:   crtol = tao->crtol;
412:   PetscCall(PetscOptionsReal("-tao_catol", "Stop if constraints violations within", "TaoSetConstraintTolerances", tao->catol, &catol, NULL));
413:   PetscCall(PetscOptionsReal("-tao_crtol", "Stop if relative constraint violations within", "TaoSetConstraintTolerances", tao->crtol, &crtol, NULL));
414:   PetscCall(TaoSetConstraintTolerances(tao, catol, crtol));

416:   gatol = tao->gatol;
417:   grtol = tao->grtol;
418:   gttol = tao->gttol;
419:   PetscCall(PetscOptionsReal("-tao_gatol", "Stop if norm of gradient less than", "TaoSetTolerances", tao->gatol, &gatol, NULL));
420:   PetscCall(PetscOptionsReal("-tao_grtol", "Stop if norm of gradient divided by the function value is less than", "TaoSetTolerances", tao->grtol, &grtol, NULL));
421:   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));
422:   PetscCall(TaoSetTolerances(tao, gatol, grtol, gttol));

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

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

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

435:   PetscCall(PetscOptionsDeprecated("-tao_solution_monitor", "-tao_monitor_solution", "3.21", NULL));
436:   PetscCall(PetscOptionsDeprecated("-tao_gradient_monitor", "-tao_monitor_gradient", "3.21", NULL));
437:   PetscCall(PetscOptionsDeprecated("-tao_stepdirection_monitor", "-tao_monitor_step", "3.21", NULL));
438:   PetscCall(PetscOptionsDeprecated("-tao_residual_monitor", "-tao_monitor_residual", "3.21", NULL));
439:   PetscCall(PetscOptionsDeprecated("-tao_smonitor", "-tao_monitor_short", "3.21", NULL));
440:   PetscCall(PetscOptionsDeprecated("-tao_cmonitor", "-tao_monitor_constraint_norm", "3.21", NULL));
441:   PetscCall(PetscOptionsDeprecated("-tao_gmonitor", "-tao_monitor_globalization", "3.21", NULL));
442:   PetscCall(PetscOptionsDeprecated("-tao_draw_solution", "-tao_monitor_solution_draw", "3.21", NULL));
443:   PetscCall(PetscOptionsDeprecated("-tao_draw_gradient", "-tao_monitor_gradient_draw", "3.21", NULL));
444:   PetscCall(PetscOptionsDeprecated("-tao_draw_step", "-tao_monitor_step_draw", "3.21", NULL));

446:   PetscCall(PetscOptionsString("-tao_monitor_solution", "View solution vector after each iteration", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
447:   if (flg) {
448:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
449:     PetscCall(TaoMonitorSet(tao, TaoMonitorSolution, monviewer, (PetscCtxDestroyFn *)PetscViewerDestroy));
450:   }

452:   PetscCall(PetscOptionsBool("-tao_converged_reason", "Print reason for Tao converged", "TaoSolve", tao->printreason, &tao->printreason, NULL));
453:   PetscCall(PetscOptionsString("-tao_monitor_gradient", "View gradient vector for each iteration", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
454:   if (flg) {
455:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
456:     PetscCall(TaoMonitorSet(tao, TaoMonitorGradient, monviewer, (PetscCtxDestroyFn *)PetscViewerDestroy));
457:   }

459:   PetscCall(PetscOptionsString("-tao_monitor_step", "View step vector after each iteration", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
460:   if (flg) {
461:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
462:     PetscCall(TaoMonitorSet(tao, TaoMonitorStep, monviewer, (PetscCtxDestroyFn *)PetscViewerDestroy));
463:   }

465:   PetscCall(PetscOptionsString("-tao_monitor_residual", "View least-squares residual vector after each iteration", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
466:   if (flg) {
467:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
468:     PetscCall(TaoMonitorSet(tao, TaoMonitorResidual, monviewer, (PetscCtxDestroyFn *)PetscViewerDestroy));
469:   }

471:   PetscCall(PetscOptionsString("-tao_monitor", "Use the default convergence monitor", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
472:   if (flg) {
473:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
474:     PetscCall(TaoMonitorSet(tao, TaoMonitorDefault, monviewer, (PetscCtxDestroyFn *)PetscViewerDestroy));
475:   }

477:   PetscCall(PetscOptionsString("-tao_monitor_globalization", "Use the convergence monitor with extra globalization info", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
478:   if (flg) {
479:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
480:     PetscCall(TaoMonitorSet(tao, TaoMonitorGlobalization, monviewer, (PetscCtxDestroyFn *)PetscViewerDestroy));
481:   }

483:   PetscCall(PetscOptionsString("-tao_monitor_short", "Use the short convergence monitor", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
484:   if (flg) {
485:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
486:     PetscCall(TaoMonitorSet(tao, TaoMonitorDefaultShort, monviewer, (PetscCtxDestroyFn *)PetscViewerDestroy));
487:   }

489:   PetscCall(PetscOptionsString("-tao_monitor_constraint_norm", "Use the default convergence monitor with constraint norm", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
490:   if (flg) {
491:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
492:     PetscCall(TaoMonitorSet(tao, TaoMonitorConstraintNorm, monviewer, (PetscCtxDestroyFn *)PetscViewerDestroy));
493:   }

495:   flg = PETSC_FALSE;
496:   PetscCall(PetscOptionsDeprecated("-tao_cancelmonitors", "-tao_monitor_cancel", "3.21", NULL));
497:   PetscCall(PetscOptionsBool("-tao_monitor_cancel", "cancel all monitors and call any registered destroy routines", "TaoMonitorCancel", flg, &flg, NULL));
498:   if (flg) PetscCall(TaoMonitorCancel(tao));

500:   flg = PETSC_FALSE;
501:   PetscCall(PetscOptionsBool("-tao_monitor_solution_draw", "Plot solution vector at each iteration", "TaoMonitorSet", flg, &flg, NULL));
502:   if (flg) {
503:     TaoMonitorDrawCtx drawctx;
504:     PetscInt          howoften = 1;
505:     PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &drawctx));
506:     PetscCall(TaoMonitorSet(tao, TaoMonitorSolutionDraw, drawctx, (PetscCtxDestroyFn *)TaoMonitorDrawCtxDestroy));
507:   }

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

513:   flg = PETSC_FALSE;
514:   PetscCall(PetscOptionsBool("-tao_monitor_gradient_draw", "plots gradient at each iteration", "TaoMonitorSet", flg, &flg, NULL));
515:   if (flg) {
516:     TaoMonitorDrawCtx drawctx;
517:     PetscInt          howoften = 1;
518:     PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &drawctx));
519:     PetscCall(TaoMonitorSet(tao, TaoMonitorGradientDraw, drawctx, (PetscCtxDestroyFn *)TaoMonitorDrawCtxDestroy));
520:   }
521:   flg = PETSC_FALSE;
522:   PetscCall(PetscOptionsBool("-tao_fd_gradient", "compute gradient using finite differences", "TaoDefaultComputeGradient", flg, &flg, NULL));
523:   if (flg) PetscCall(TaoSetGradient(tao, NULL, TaoDefaultComputeGradient, NULL));
524:   flg = PETSC_FALSE;
525:   PetscCall(PetscOptionsBool("-tao_fd_hessian", "compute Hessian using finite differences", "TaoDefaultComputeHessian", flg, &flg, NULL));
526:   if (flg) {
527:     Mat H;

529:     PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &H));
530:     PetscCall(MatSetType(H, MATAIJ));
531:     PetscCall(TaoSetHessian(tao, H, H, TaoDefaultComputeHessian, NULL));
532:     PetscCall(MatDestroy(&H));
533:   }
534:   flg = PETSC_FALSE;
535:   PetscCall(PetscOptionsBool("-tao_mf_hessian", "compute matrix-free Hessian using finite differences", "TaoDefaultComputeHessianMFFD", flg, &flg, NULL));
536:   if (flg) {
537:     Mat H;

539:     PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &H));
540:     PetscCall(TaoSetHessian(tao, H, H, TaoDefaultComputeHessianMFFD, NULL));
541:     PetscCall(MatDestroy(&H));
542:   }
543:   PetscCall(PetscOptionsBool("-tao_recycle_history", "enable recycling/re-using information from the previous TaoSolve() call for some algorithms", "TaoSetRecycleHistory", flg, &flg, &found));
544:   if (found) PetscCall(TaoSetRecycleHistory(tao, flg));
545:   PetscCall(PetscOptionsEnum("-tao_subset_type", "subset type", "", TaoSubSetTypes, (PetscEnum)tao->subset_type, (PetscEnum *)&tao->subset_type, NULL));

547:   if (tao->ksp) {
548:     PetscCall(PetscOptionsBool("-tao_ksp_ew", "Use Eisentat-Walker linear system convergence test", "TaoKSPSetUseEW", tao->ksp_ewconv, &tao->ksp_ewconv, NULL));
549:     PetscCall(TaoKSPSetUseEW(tao, tao->ksp_ewconv));
550:   }

552:   PetscTryTypeMethod(tao, setfromoptions, PetscOptionsObject);

554:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
555:   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)tao, PetscOptionsObject));
556:   PetscOptionsEnd();

558:   if (tao->linesearch) PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
559:   PetscFunctionReturn(PETSC_SUCCESS);
560: }

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

565:   Collective

567:   Input Parameters:
568: + A    - the  `Tao` context
569: . obj  - Optional object that provides the prefix for the options database
570: - name - command line option

572:   Level: intermediate

574: .seealso: [](ch_tao), `Tao`, `TaoView`, `PetscObjectViewFromOptions()`, `TaoCreate()`
575: @*/
576: PetscErrorCode TaoViewFromOptions(Tao A, PetscObject obj, const char name[])
577: {
578:   PetscFunctionBegin;
580:   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
581:   PetscFunctionReturn(PETSC_SUCCESS);
582: }

584: /*@
585:   TaoView - Prints information about the `Tao` object

587:   Collective

589:   Input Parameters:
590: + tao    - the `Tao` context
591: - viewer - visualization context

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

596:   Level: beginner

598:   Notes:
599:   The available visualization contexts include
600: +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
601: -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
602:   output where only the first processor opens
603:   the file.  All other processors send their
604:   data to the first processor to print.

606: .seealso: [](ch_tao), `Tao`, `PetscViewerASCIIOpen()`
607: @*/
608: PetscErrorCode TaoView(Tao tao, PetscViewer viewer)
609: {
610:   PetscBool isascii, isstring;
611:   TaoType   type;

613:   PetscFunctionBegin;
615:   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(((PetscObject)tao)->comm, &viewer));
617:   PetscCheckSameComm(tao, 1, viewer, 2);

619:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
620:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
621:   if (isascii) {
622:     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)tao, viewer));

624:     PetscCall(PetscViewerASCIIPushTab(viewer));
625:     PetscTryTypeMethod(tao, view, viewer);
626:     if (tao->linesearch) PetscCall(TaoLineSearchView(tao->linesearch, viewer));
627:     if (tao->ksp) {
628:       PetscCall(KSPView(tao->ksp, viewer));
629:       PetscCall(PetscViewerASCIIPrintf(viewer, "total KSP iterations: %" PetscInt_FMT "\n", tao->ksp_tot_its));
630:     }

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

634:     PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: gatol=%g,", (double)tao->gatol));
635:     PetscCall(PetscViewerASCIIPrintf(viewer, " grtol=%g,", (double)tao->grtol));
636:     PetscCall(PetscViewerASCIIPrintf(viewer, " steptol=%g,", (double)tao->steptol));
637:     PetscCall(PetscViewerASCIIPrintf(viewer, " gttol=%g\n", (double)tao->gttol));
638:     PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Function/Gradient:=%g\n", (double)tao->residual));

640:     if (tao->constrained) {
641:       PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances:"));
642:       PetscCall(PetscViewerASCIIPrintf(viewer, " catol=%g,", (double)tao->catol));
643:       PetscCall(PetscViewerASCIIPrintf(viewer, " crtol=%g\n", (double)tao->crtol));
644:       PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Constraints:=%g\n", (double)tao->cnorm));
645:     }

647:     if (tao->trust < tao->steptol) {
648:       PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: steptol=%g\n", (double)tao->steptol));
649:       PetscCall(PetscViewerASCIIPrintf(viewer, "Final trust region radius:=%g\n", (double)tao->trust));
650:     }

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

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

658:     if (tao->nfuncs > 0) {
659:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function evaluations=%" PetscInt_FMT ",", tao->nfuncs));
660:       if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: unlimited)\n"));
661:       else PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: %" PetscInt_FMT ")\n", tao->max_funcs));
662:     }
663:     if (tao->ngrads > 0) {
664:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of gradient evaluations=%" PetscInt_FMT ",", tao->ngrads));
665:       if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: unlimited)\n"));
666:       else PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: %" PetscInt_FMT ")\n", tao->max_funcs));
667:     }
668:     if (tao->nfuncgrads > 0) {
669:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function/gradient evaluations=%" PetscInt_FMT ",", tao->nfuncgrads));
670:       if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, "    (max: unlimited)\n"));
671:       else PetscCall(PetscViewerASCIIPrintf(viewer, "    (max: %" PetscInt_FMT ")\n", tao->max_funcs));
672:     }
673:     if (tao->nhess > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Hessian evaluations=%" PetscInt_FMT "\n", tao->nhess));
674:     if (tao->nconstraints > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of constraint function evaluations=%" PetscInt_FMT "\n", tao->nconstraints));
675:     if (tao->njac > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Jacobian evaluations=%" PetscInt_FMT "\n", tao->njac));

677:     if (tao->reason > 0) {
678:       PetscCall(PetscViewerASCIIPrintf(viewer, "Solution converged: "));
679:       switch (tao->reason) {
680:       case TAO_CONVERGED_GATOL:
681:         PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)|| <= gatol\n"));
682:         break;
683:       case TAO_CONVERGED_GRTOL:
684:         PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)||/|f(X)| <= grtol\n"));
685:         break;
686:       case TAO_CONVERGED_GTTOL:
687:         PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)||/||g(X0)|| <= gttol\n"));
688:         break;
689:       case TAO_CONVERGED_STEPTOL:
690:         PetscCall(PetscViewerASCIIPrintf(viewer, " Steptol -- step size small\n"));
691:         break;
692:       case TAO_CONVERGED_MINF:
693:         PetscCall(PetscViewerASCIIPrintf(viewer, " Minf --  f < fmin\n"));
694:         break;
695:       case TAO_CONVERGED_USER:
696:         PetscCall(PetscViewerASCIIPrintf(viewer, " User Terminated\n"));
697:         break;
698:       default:
699:         PetscCall(PetscViewerASCIIPrintf(viewer, " %d\n", tao->reason));
700:         break;
701:       }
702:     } else if (tao->reason == TAO_CONTINUE_ITERATING) {
703:       PetscCall(PetscViewerASCIIPrintf(viewer, "Solver never run\n"));
704:     } else {
705:       PetscCall(PetscViewerASCIIPrintf(viewer, "Solver failed: "));
706:       switch (tao->reason) {
707:       case TAO_DIVERGED_MAXITS:
708:         PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum Iterations\n"));
709:         break;
710:       case TAO_DIVERGED_NAN:
711:         PetscCall(PetscViewerASCIIPrintf(viewer, " NAN or Inf encountered\n"));
712:         break;
713:       case TAO_DIVERGED_MAXFCN:
714:         PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum Function Evaluations\n"));
715:         break;
716:       case TAO_DIVERGED_LS_FAILURE:
717:         PetscCall(PetscViewerASCIIPrintf(viewer, " Line Search Failure\n"));
718:         break;
719:       case TAO_DIVERGED_TR_REDUCTION:
720:         PetscCall(PetscViewerASCIIPrintf(viewer, " Trust Region too small\n"));
721:         break;
722:       case TAO_DIVERGED_USER:
723:         PetscCall(PetscViewerASCIIPrintf(viewer, " User Terminated\n"));
724:         break;
725:       default:
726:         PetscCall(PetscViewerASCIIPrintf(viewer, " %d\n", tao->reason));
727:         break;
728:       }
729:     }
730:     PetscCall(PetscViewerASCIIPopTab(viewer));
731:   } else if (isstring) {
732:     PetscCall(TaoGetType(tao, &type));
733:     PetscCall(PetscViewerStringSPrintf(viewer, " %-3.3s", type));
734:   }
735:   PetscFunctionReturn(PETSC_SUCCESS);
736: }

738: /*@
739:   TaoSetRecycleHistory - Sets the boolean flag to enable/disable re-using
740:   iterate information from the previous `TaoSolve()`. This feature is disabled by
741:   default.

743:   Logically Collective

745:   Input Parameters:
746: + tao     - the `Tao` context
747: - recycle - boolean flag

749:   Options Database Key:
750: . -tao_recycle_history <true,false> - reuse the history

752:   Level: intermediate

754:   Notes:
755:   For conjugate gradient methods (`TAOBNCG`), this re-uses the latest search direction
756:   from the previous `TaoSolve()` call when computing the first search direction in a
757:   new solution. By default, CG methods set the first search direction to the
758:   negative gradient.

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

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

767: .seealso: [](ch_tao), `Tao`, `TaoGetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
768: @*/
769: PetscErrorCode TaoSetRecycleHistory(Tao tao, PetscBool recycle)
770: {
771:   PetscFunctionBegin;
774:   tao->recycle = recycle;
775:   PetscFunctionReturn(PETSC_SUCCESS);
776: }

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

782:   Logically Collective

784:   Input Parameter:
785: . tao - the `Tao` context

787:   Output Parameter:
788: . recycle - boolean flag

790:   Level: intermediate

792: .seealso: [](ch_tao), `Tao`, `TaoSetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
793: @*/
794: PetscErrorCode TaoGetRecycleHistory(Tao tao, PetscBool *recycle)
795: {
796:   PetscFunctionBegin;
798:   PetscAssertPointer(recycle, 2);
799:   *recycle = tao->recycle;
800:   PetscFunctionReturn(PETSC_SUCCESS);
801: }

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

806:   Logically Collective

808:   Input Parameters:
809: + tao   - the `Tao` context
810: . gatol - stop if norm of gradient is less than this
811: . grtol - stop if relative norm of gradient is less than this
812: - gttol - stop if norm of gradient is reduced by this factor

814:   Options Database Keys:
815: + -tao_gatol <gatol> - Sets gatol
816: . -tao_grtol <grtol> - Sets grtol
817: - -tao_gttol <gttol> - Sets gttol

819:   Stopping Criteria\:
820: .vb
821:   ||g(X)||                            <= gatol
822:   ||g(X)|| / |f(X)|                   <= grtol
823:   ||g(X)|| / ||g(X0)||                <= gttol
824: .ve

826:   Level: beginner

828:   Notes:
829:   Use `PETSC_CURRENT` to leave one or more tolerances unchanged.

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

833:   Fortran Note:
834:   Use `PETSC_CURRENT_REAL` or `PETSC_DETERMINE_REAL`

836: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`
837: @*/
838: PetscErrorCode TaoSetTolerances(Tao tao, PetscReal gatol, PetscReal grtol, PetscReal gttol)
839: {
840:   PetscFunctionBegin;

846:   if (gatol == (PetscReal)PETSC_DETERMINE) {
847:     tao->gatol = tao->default_gatol;
848:   } else if (gatol != (PetscReal)PETSC_CURRENT) {
849:     PetscCheck(gatol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative gatol not allowed");
850:     tao->gatol = gatol;
851:   }

853:   if (grtol == (PetscReal)PETSC_DETERMINE) {
854:     tao->grtol = tao->default_grtol;
855:   } else if (grtol != (PetscReal)PETSC_CURRENT) {
856:     PetscCheck(grtol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative grtol not allowed");
857:     tao->grtol = grtol;
858:   }

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

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

872:   Logically Collective

874:   Input Parameters:
875: + tao   - the `Tao` context
876: . catol - absolute constraint tolerance, constraint norm must be less than `catol` for used for `gatol` convergence criteria
877: - crtol - relative constraint tolerance, constraint norm must be less than `crtol` for used for `gatol`, `gttol` convergence criteria

879:   Options Database Keys:
880: + -tao_catol <catol> - Sets catol
881: - -tao_crtol <crtol> - Sets crtol

883:   Level: intermediate

885:   Notes:
886:   Use `PETSC_CURRENT` to leave one or tolerance unchanged.

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

890:   Fortran Note:
891:   Use `PETSC_CURRENT_REAL` or `PETSC_DETERMINE_REAL`

893: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`, `TaoGetConstraintTolerances()`, `TaoSetTolerances()`
894: @*/
895: PetscErrorCode TaoSetConstraintTolerances(Tao tao, PetscReal catol, PetscReal crtol)
896: {
897:   PetscFunctionBegin;

902:   if (catol == (PetscReal)PETSC_DETERMINE) {
903:     tao->catol = tao->default_catol;
904:   } else if (catol != (PetscReal)PETSC_CURRENT) {
905:     PetscCheck(catol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative catol not allowed");
906:     tao->catol = catol;
907:   }

909:   if (crtol == (PetscReal)PETSC_DETERMINE) {
910:     tao->crtol = tao->default_crtol;
911:   } else if (crtol != (PetscReal)PETSC_CURRENT) {
912:     PetscCheck(crtol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative crtol not allowed");
913:     tao->crtol = crtol;
914:   }
915:   PetscFunctionReturn(PETSC_SUCCESS);
916: }

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

921:   Not Collective

923:   Input Parameter:
924: . tao - the `Tao` context

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

930:   Level: intermediate

932: .seealso: [](ch_tao), `Tao`, `TaoConvergedReasons`,`TaoGetTolerances()`, `TaoSetTolerances()`, `TaoSetConstraintTolerances()`
933: @*/
934: PetscErrorCode TaoGetConstraintTolerances(Tao tao, PetscReal *catol, PetscReal *crtol)
935: {
936:   PetscFunctionBegin;
938:   if (catol) *catol = tao->catol;
939:   if (crtol) *crtol = tao->crtol;
940:   PetscFunctionReturn(PETSC_SUCCESS);
941: }

943: /*@
944:   TaoSetFunctionLowerBound - Sets a bound on the solution objective value.
945:   When an approximate solution with an objective value below this number
946:   has been found, the solver will terminate.

948:   Logically Collective

950:   Input Parameters:
951: + tao  - the Tao solver context
952: - fmin - the tolerance

954:   Options Database Key:
955: . -tao_fmin <fmin> - sets the minimum function value

957:   Level: intermediate

959: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetTolerances()`
960: @*/
961: PetscErrorCode TaoSetFunctionLowerBound(Tao tao, PetscReal fmin)
962: {
963:   PetscFunctionBegin;
966:   tao->fmin = fmin;
967:   PetscFunctionReturn(PETSC_SUCCESS);
968: }

970: /*@
971:   TaoGetFunctionLowerBound - Gets the bound on the solution objective value.
972:   When an approximate solution with an objective value below this number
973:   has been found, the solver will terminate.

975:   Not Collective

977:   Input Parameter:
978: . tao - the `Tao` solver context

980:   Output Parameter:
981: . fmin - the minimum function value

983:   Level: intermediate

985: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetFunctionLowerBound()`
986: @*/
987: PetscErrorCode TaoGetFunctionLowerBound(Tao tao, PetscReal *fmin)
988: {
989:   PetscFunctionBegin;
991:   PetscAssertPointer(fmin, 2);
992:   *fmin = tao->fmin;
993:   PetscFunctionReturn(PETSC_SUCCESS);
994: }

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

999:   Logically Collective

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

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

1008:   Level: intermediate

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

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

1016: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoSetMaximumIterations()`
1017: @*/
1018: PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao tao, PetscInt nfcn)
1019: {
1020:   PetscFunctionBegin;
1023:   if (nfcn == PETSC_DETERMINE) {
1024:     tao->max_funcs = tao->default_max_funcs;
1025:   } else if (nfcn == PETSC_UNLIMITED || nfcn < 0) {
1026:     tao->max_funcs = PETSC_UNLIMITED;
1027:   } else {
1028:     PetscCheck(nfcn >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of function evaluations  must be positive");
1029:     tao->max_funcs = nfcn;
1030:   }
1031:   PetscFunctionReturn(PETSC_SUCCESS);
1032: }

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

1037:   Logically Collective

1039:   Input Parameter:
1040: . tao - the `Tao` solver context

1042:   Output Parameter:
1043: . nfcn - the maximum number of function evaluations

1045:   Level: intermediate

1047: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1048: @*/
1049: PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao tao, PetscInt *nfcn)
1050: {
1051:   PetscFunctionBegin;
1053:   PetscAssertPointer(nfcn, 2);
1054:   *nfcn = tao->max_funcs;
1055:   PetscFunctionReturn(PETSC_SUCCESS);
1056: }

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

1061:   Not Collective

1063:   Input Parameter:
1064: . tao - the `Tao` solver context

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

1069:   Level: intermediate

1071: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1072: @*/
1073: PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao tao, PetscInt *nfuncs)
1074: {
1075:   PetscFunctionBegin;
1077:   PetscAssertPointer(nfuncs, 2);
1078:   *nfuncs = PetscMax(tao->nfuncs, tao->nfuncgrads);
1079:   PetscFunctionReturn(PETSC_SUCCESS);
1080: }

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

1085:   Logically Collective

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

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

1094:   Level: intermediate

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

1099:   Developer Note:
1100:   DeprAlso accepts the deprecated negative values to indicate no limit

1102: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoSetMaximumFunctionEvaluations()`
1103: @*/
1104: PetscErrorCode TaoSetMaximumIterations(Tao tao, PetscInt maxits)
1105: {
1106:   PetscFunctionBegin;
1109:   if (maxits == PETSC_DETERMINE) {
1110:     tao->max_it = tao->default_max_it;
1111:   } else if (maxits == PETSC_UNLIMITED) {
1112:     tao->max_it = PETSC_INT_MAX;
1113:   } else {
1114:     PetscCheck(maxits > 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations must be positive");
1115:     tao->max_it = maxits;
1116:   }
1117:   PetscFunctionReturn(PETSC_SUCCESS);
1118: }

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

1123:   Not Collective

1125:   Input Parameter:
1126: . tao - the `Tao` solver context

1128:   Output Parameter:
1129: . maxits - the maximum number of iterates

1131:   Level: intermediate

1133: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumIterations()`, `TaoGetMaximumFunctionEvaluations()`
1134: @*/
1135: PetscErrorCode TaoGetMaximumIterations(Tao tao, PetscInt *maxits)
1136: {
1137:   PetscFunctionBegin;
1139:   PetscAssertPointer(maxits, 2);
1140:   *maxits = tao->max_it;
1141:   PetscFunctionReturn(PETSC_SUCCESS);
1142: }

1144: /*@
1145:   TaoSetInitialTrustRegionRadius - Sets the initial trust region radius.

1147:   Logically Collective

1149:   Input Parameters:
1150: + tao    - a `Tao` optimization solver
1151: - radius - the trust region radius

1153:   Options Database Key:
1154: . -tao_trust0 <t0> - sets initial trust region radius

1156:   Level: intermediate

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

1161: .seealso: [](ch_tao), `Tao`, `TaoGetTrustRegionRadius()`, `TaoSetTrustRegionTolerance()`, `TAONTR`
1162: @*/
1163: PetscErrorCode TaoSetInitialTrustRegionRadius(Tao tao, PetscReal radius)
1164: {
1165:   PetscFunctionBegin;
1168:   if (radius == PETSC_DETERMINE) {
1169:     tao->trust0 = tao->default_trust0;
1170:   } else {
1171:     PetscCheck(radius > 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Radius must be positive");
1172:     tao->trust0 = radius;
1173:   }
1174:   PetscFunctionReturn(PETSC_SUCCESS);
1175: }

1177: /*@
1178:   TaoGetInitialTrustRegionRadius - Gets the initial trust region radius.

1180:   Not Collective

1182:   Input Parameter:
1183: . tao - a `Tao` optimization solver

1185:   Output Parameter:
1186: . radius - the trust region radius

1188:   Level: intermediate

1190: .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetCurrentTrustRegionRadius()`, `TAONTR`
1191: @*/
1192: PetscErrorCode TaoGetInitialTrustRegionRadius(Tao tao, PetscReal *radius)
1193: {
1194:   PetscFunctionBegin;
1196:   PetscAssertPointer(radius, 2);
1197:   *radius = tao->trust0;
1198:   PetscFunctionReturn(PETSC_SUCCESS);
1199: }

1201: /*@
1202:   TaoGetCurrentTrustRegionRadius - Gets the current trust region radius.

1204:   Not Collective

1206:   Input Parameter:
1207: . tao - a `Tao` optimization solver

1209:   Output Parameter:
1210: . radius - the trust region radius

1212:   Level: intermediate

1214: .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetInitialTrustRegionRadius()`, `TAONTR`
1215: @*/
1216: PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao tao, PetscReal *radius)
1217: {
1218:   PetscFunctionBegin;
1220:   PetscAssertPointer(radius, 2);
1221:   *radius = tao->trust;
1222:   PetscFunctionReturn(PETSC_SUCCESS);
1223: }

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

1228:   Not Collective

1230:   Input Parameter:
1231: . tao - the `Tao` context

1233:   Output Parameters:
1234: + gatol - stop if norm of gradient is less than this
1235: . grtol - stop if relative norm of gradient is less than this
1236: - gttol - stop if norm of gradient is reduced by a this factor

1238:   Level: intermediate

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

1243: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`
1244: @*/
1245: PetscErrorCode TaoGetTolerances(Tao tao, PetscReal *gatol, PetscReal *grtol, PetscReal *gttol)
1246: {
1247:   PetscFunctionBegin;
1249:   if (gatol) *gatol = tao->gatol;
1250:   if (grtol) *grtol = tao->grtol;
1251:   if (gttol) *gttol = tao->gttol;
1252:   PetscFunctionReturn(PETSC_SUCCESS);
1253: }

1255: /*@
1256:   TaoGetKSP - Gets the linear solver used by the optimization solver.

1258:   Not Collective

1260:   Input Parameter:
1261: . tao - the `Tao` solver

1263:   Output Parameter:
1264: . ksp - the `KSP` linear solver used in the optimization solver

1266:   Level: intermediate

1268: .seealso: [](ch_tao), `Tao`, `KSP`
1269: @*/
1270: PetscErrorCode TaoGetKSP(Tao tao, KSP *ksp)
1271: {
1272:   PetscFunctionBegin;
1274:   PetscAssertPointer(ksp, 2);
1275:   *ksp = tao->ksp;
1276:   PetscFunctionReturn(PETSC_SUCCESS);
1277: }

1279: /*@
1280:   TaoGetLinearSolveIterations - Gets the total number of linear iterations
1281:   used by the `Tao` solver

1283:   Not Collective

1285:   Input Parameter:
1286: . tao - the `Tao` context

1288:   Output Parameter:
1289: . lits - number of linear iterations

1291:   Level: intermediate

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

1296: .seealso: [](ch_tao), `Tao`, `TaoGetKSP()`
1297: @*/
1298: PetscErrorCode TaoGetLinearSolveIterations(Tao tao, PetscInt *lits)
1299: {
1300:   PetscFunctionBegin;
1302:   PetscAssertPointer(lits, 2);
1303:   *lits = tao->ksp_tot_its;
1304:   PetscFunctionReturn(PETSC_SUCCESS);
1305: }

1307: /*@
1308:   TaoGetLineSearch - Gets the line search used by the optimization solver.

1310:   Not Collective

1312:   Input Parameter:
1313: . tao - the `Tao` solver

1315:   Output Parameter:
1316: . ls - the line search used in the optimization solver

1318:   Level: intermediate

1320: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`
1321: @*/
1322: PetscErrorCode TaoGetLineSearch(Tao tao, TaoLineSearch *ls)
1323: {
1324:   PetscFunctionBegin;
1326:   PetscAssertPointer(ls, 2);
1327:   *ls = tao->linesearch;
1328:   PetscFunctionReturn(PETSC_SUCCESS);
1329: }

1331: /*@
1332:   TaoAddLineSearchCounts - Adds the number of function evaluations spent
1333:   in the line search to the running total.

1335:   Input Parameters:
1336: . tao - the `Tao` solver

1338:   Level: developer

1340: .seealso: [](ch_tao), `Tao`, `TaoGetLineSearch()`, `TaoLineSearchApply()`
1341: @*/
1342: PetscErrorCode TaoAddLineSearchCounts(Tao tao)
1343: {
1344:   PetscBool flg;
1345:   PetscInt  nfeval, ngeval, nfgeval;

1347:   PetscFunctionBegin;
1349:   if (tao->linesearch) {
1350:     PetscCall(TaoLineSearchIsUsingTaoRoutines(tao->linesearch, &flg));
1351:     if (!flg) {
1352:       PetscCall(TaoLineSearchGetNumberFunctionEvaluations(tao->linesearch, &nfeval, &ngeval, &nfgeval));
1353:       tao->nfuncs += nfeval;
1354:       tao->ngrads += ngeval;
1355:       tao->nfuncgrads += nfgeval;
1356:     }
1357:   }
1358:   PetscFunctionReturn(PETSC_SUCCESS);
1359: }

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

1364:   Not Collective

1366:   Input Parameter:
1367: . tao - the `Tao` context

1369:   Output Parameter:
1370: . X - the current solution

1372:   Level: intermediate

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

1377: .seealso: [](ch_tao), `Tao`, `TaoSetSolution()`, `TaoSolve()`
1378: @*/
1379: PetscErrorCode TaoGetSolution(Tao tao, Vec *X)
1380: {
1381:   PetscFunctionBegin;
1383:   PetscAssertPointer(X, 2);
1384:   *X = tao->solution;
1385:   PetscFunctionReturn(PETSC_SUCCESS);
1386: }

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

1393:   Collective

1395:   Input Parameter:
1396: . tao - the `Tao` context

1398:   Level: developer

1400: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
1401: @*/
1402: PetscErrorCode TaoResetStatistics(Tao tao)
1403: {
1404:   PetscFunctionBegin;
1406:   tao->niter        = 0;
1407:   tao->nfuncs       = 0;
1408:   tao->nfuncgrads   = 0;
1409:   tao->ngrads       = 0;
1410:   tao->nhess        = 0;
1411:   tao->njac         = 0;
1412:   tao->nconstraints = 0;
1413:   tao->ksp_its      = 0;
1414:   tao->ksp_tot_its  = 0;
1415:   tao->reason       = TAO_CONTINUE_ITERATING;
1416:   tao->residual     = 0.0;
1417:   tao->cnorm        = 0.0;
1418:   tao->step         = 0.0;
1419:   tao->lsflag       = PETSC_FALSE;
1420:   if (tao->hist_reset) tao->hist_len = 0;
1421:   PetscFunctionReturn(PETSC_SUCCESS);
1422: }

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

1429:   Logically Collective

1431:   Input Parameters:
1432: + tao  - The `Tao` solver
1433: . func - The function
1434: - ctx  - The update function context

1436:   Calling sequence of `func`:
1437: + tao - The optimizer context
1438: . it  - The current iteration index
1439: - ctx - The update context

1441:   Level: advanced

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

1447: .seealso: [](ch_tao), `Tao`, `TaoSolve()`
1448: @*/
1449: PetscErrorCode TaoSetUpdate(Tao tao, PetscErrorCode (*func)(Tao tao, PetscInt it, void *ctx), void *ctx)
1450: {
1451:   PetscFunctionBegin;
1453:   tao->ops->update = func;
1454:   tao->user_update = ctx;
1455:   PetscFunctionReturn(PETSC_SUCCESS);
1456: }

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

1463:   Logically Collective

1465:   Input Parameters:
1466: + tao  - the `Tao` object
1467: . conv - the routine to test for convergence
1468: - ctx  - [optional] context for private data for the convergence routine
1469:         (may be `NULL`)

1471:   Calling sequence of `conv`:
1472: + tao - the `Tao` object
1473: - ctx - [optional] convergence context

1475:   Level: advanced

1477:   Note:
1478:   The new convergence testing routine should call `TaoSetConvergedReason()`.

1480: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetConvergedReason()`, `TaoGetSolutionStatus()`, `TaoGetTolerances()`, `TaoMonitorSet()`
1481: @*/
1482: PetscErrorCode TaoSetConvergenceTest(Tao tao, PetscErrorCode (*conv)(Tao, void *), void *ctx)
1483: {
1484:   PetscFunctionBegin;
1486:   tao->ops->convergencetest = conv;
1487:   tao->cnvP                 = ctx;
1488:   PetscFunctionReturn(PETSC_SUCCESS);
1489: }

1491: /*@C
1492:   TaoMonitorSet - Sets an additional function that is to be used at every
1493:   iteration of the solver to display the iteration's
1494:   progress.

1496:   Logically Collective

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

1504:   Calling sequence of `func`:
1505: + tao - the `Tao` solver context
1506: - ctx - [optional] monitoring context

1508:   Level: intermediate

1510:   Notes:
1511:   See `TaoSetFromOptions()` for a monitoring options.

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

1517:   Fortran Notes:
1518:   Only one monitor function may be set

1520: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoMonitorDefault()`, `TaoMonitorCancel()`, `TaoSetDestroyRoutine()`, `TaoView()`, `PetscCtxDestroyFn`
1521: @*/
1522: PetscErrorCode TaoMonitorSet(Tao tao, PetscErrorCode (*func)(Tao, void *), void *ctx, PetscCtxDestroyFn *dest)
1523: {
1524:   PetscInt  i;
1525:   PetscBool identical;

1527:   PetscFunctionBegin;
1529:   PetscCheck(tao->numbermonitors < MAXTAOMONITORS, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Cannot attach another monitor -- max=%d", MAXTAOMONITORS);

1531:   for (i = 0; i < tao->numbermonitors; i++) {
1532:     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))func, ctx, dest, (PetscErrorCode (*)(void))tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i], &identical));
1533:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
1534:   }
1535:   tao->monitor[tao->numbermonitors]        = func;
1536:   tao->monitorcontext[tao->numbermonitors] = ctx;
1537:   tao->monitordestroy[tao->numbermonitors] = dest;
1538:   ++tao->numbermonitors;
1539:   PetscFunctionReturn(PETSC_SUCCESS);
1540: }

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

1545:   Logically Collective

1547:   Input Parameter:
1548: . tao - the `Tao` solver context

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

1555:   Level: advanced

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

1560: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1561: @*/
1562: PetscErrorCode TaoMonitorCancel(Tao tao)
1563: {
1564:   PetscInt i;

1566:   PetscFunctionBegin;
1568:   for (i = 0; i < tao->numbermonitors; i++) {
1569:     if (tao->monitordestroy[i]) PetscCall((*tao->monitordestroy[i])(&tao->monitorcontext[i]));
1570:   }
1571:   tao->numbermonitors = 0;
1572:   PetscFunctionReturn(PETSC_SUCCESS);
1573: }

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

1578:   Collective

1580:   Input Parameters:
1581: + tao - the `Tao` context
1582: - ctx - `PetscViewer` context or `NULL`

1584:   Options Database Key:
1585: . -tao_monitor - turn on default monitoring

1587:   Level: advanced

1589:   Note:
1590:   This monitor prints the function value and gradient
1591:   norm at each iteration.

1593: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1594: @*/
1595: PetscErrorCode TaoMonitorDefault(Tao tao, void *ctx)
1596: {
1597:   PetscInt    its, tabs;
1598:   PetscReal   fct, gnorm;
1599:   PetscViewer viewer = (PetscViewer)ctx;

1601:   PetscFunctionBegin;
1604:   its   = tao->niter;
1605:   fct   = tao->fc;
1606:   gnorm = tao->residual;
1607:   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1608:   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1609:   if (its == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
1610:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Iteration information for %s solve.\n", ((PetscObject)tao)->prefix));
1611:     tao->header_printed = PETSC_TRUE;
1612:   }
1613:   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " TAO,", its));
1614:   PetscCall(PetscViewerASCIIPrintf(viewer, "  Function value: %g,", (double)fct));
1615:   if (gnorm >= PETSC_INFINITY) {
1616:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: Inf \n"));
1617:   } else {
1618:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g \n", (double)gnorm));
1619:   }
1620:   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1621:   PetscFunctionReturn(PETSC_SUCCESS);
1622: }

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

1627:   Collective

1629:   Input Parameters:
1630: + tao - the `Tao` context
1631: - ctx - `PetscViewer` context or `NULL`

1633:   Options Database Key:
1634: . -tao_monitor_globalization - turn on monitoring with globalization information

1636:   Level: advanced

1638:   Note:
1639:   This monitor prints the function value and gradient norm at each
1640:   iteration, as well as the step size and trust radius. Note that the
1641:   step size and trust radius may be the same for some algorithms.

1643: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1644: @*/
1645: PetscErrorCode TaoMonitorGlobalization(Tao tao, void *ctx)
1646: {
1647:   PetscInt    its, tabs;
1648:   PetscReal   fct, gnorm, stp, tr;
1649:   PetscViewer viewer = (PetscViewer)ctx;

1651:   PetscFunctionBegin;
1654:   its   = tao->niter;
1655:   fct   = tao->fc;
1656:   gnorm = tao->residual;
1657:   stp   = tao->step;
1658:   tr    = tao->trust;
1659:   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1660:   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1661:   if (its == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
1662:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Iteration information for %s solve.\n", ((PetscObject)tao)->prefix));
1663:     tao->header_printed = PETSC_TRUE;
1664:   }
1665:   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " TAO,", its));
1666:   PetscCall(PetscViewerASCIIPrintf(viewer, "  Function value: %g,", (double)fct));
1667:   if (gnorm >= PETSC_INFINITY) {
1668:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: Inf,"));
1669:   } else {
1670:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g,", (double)gnorm));
1671:   }
1672:   PetscCall(PetscViewerASCIIPrintf(viewer, "  Step: %g,  Trust: %g\n", (double)stp, (double)tr));
1673:   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1674:   PetscFunctionReturn(PETSC_SUCCESS);
1675: }

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

1680:   Collective

1682:   Input Parameters:
1683: + tao - the `Tao` context
1684: - ctx - `PetscViewer` context of type `PETSCVIEWERASCII`

1686:   Options Database Key:
1687: . -tao_monitor_short - turn on default short monitoring

1689:   Level: advanced

1691:   Note:
1692:   Same as `TaoMonitorDefault()` except
1693:   it prints fewer digits of the residual as the residual gets smaller.
1694:   This is because the later digits are meaningless and are often
1695:   different on different machines; by using this routine different
1696:   machines will usually generate the same output.

1698: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1699: @*/
1700: PetscErrorCode TaoMonitorDefaultShort(Tao tao, void *ctx)
1701: {
1702:   PetscInt    its, tabs;
1703:   PetscReal   fct, gnorm;
1704:   PetscViewer viewer = (PetscViewer)ctx;

1706:   PetscFunctionBegin;
1709:   its   = tao->niter;
1710:   fct   = tao->fc;
1711:   gnorm = tao->residual;
1712:   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1713:   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1714:   PetscCall(PetscViewerASCIIPrintf(viewer, "iter = %3" PetscInt_FMT ",", its));
1715:   PetscCall(PetscViewerASCIIPrintf(viewer, " Function value %g,", (double)fct));
1716:   if (gnorm >= PETSC_INFINITY) {
1717:     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: Inf \n"));
1718:   } else if (gnorm > 1.e-6) {
1719:     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: %g \n", (double)gnorm));
1720:   } else if (gnorm > 1.e-11) {
1721:     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: < 1.0e-6 \n"));
1722:   } else {
1723:     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: < 1.0e-11 \n"));
1724:   }
1725:   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1726:   PetscFunctionReturn(PETSC_SUCCESS);
1727: }

1729: /*@
1730:   TaoMonitorConstraintNorm - same as `TaoMonitorDefault()` except
1731:   it prints the norm of the constraint function.

1733:   Collective

1735:   Input Parameters:
1736: + tao - the `Tao` context
1737: - ctx - `PetscViewer` context or `NULL`

1739:   Options Database Key:
1740: . -tao_monitor_constraint_norm - monitor the constraints

1742:   Level: advanced

1744: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1745: @*/
1746: PetscErrorCode TaoMonitorConstraintNorm(Tao tao, void *ctx)
1747: {
1748:   PetscInt    its, tabs;
1749:   PetscReal   fct, gnorm;
1750:   PetscViewer viewer = (PetscViewer)ctx;

1752:   PetscFunctionBegin;
1755:   its   = tao->niter;
1756:   fct   = tao->fc;
1757:   gnorm = tao->residual;
1758:   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1759:   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1760:   PetscCall(PetscViewerASCIIPrintf(viewer, "iter = %" PetscInt_FMT ",", its));
1761:   PetscCall(PetscViewerASCIIPrintf(viewer, " Function value: %g,", (double)fct));
1762:   PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g ", (double)gnorm));
1763:   PetscCall(PetscViewerASCIIPrintf(viewer, "  Constraint: %g \n", (double)tao->cnorm));
1764:   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1765:   PetscFunctionReturn(PETSC_SUCCESS);
1766: }

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

1771:   Collective

1773:   Input Parameters:
1774: + tao - the `Tao` context
1775: - ctx - `PetscViewer` context or `NULL`

1777:   Options Database Key:
1778: . -tao_monitor_solution - view the solution

1780:   Level: advanced

1782: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1783: @*/
1784: PetscErrorCode TaoMonitorSolution(Tao tao, void *ctx)
1785: {
1786:   PetscViewer viewer = (PetscViewer)ctx;

1788:   PetscFunctionBegin;
1791:   PetscCall(VecView(tao->solution, viewer));
1792:   PetscFunctionReturn(PETSC_SUCCESS);
1793: }

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

1798:   Collective

1800:   Input Parameters:
1801: + tao - the `Tao` context
1802: - ctx - `PetscViewer` context or `NULL`

1804:   Options Database Key:
1805: . -tao_monitor_gradient - view the gradient at each iteration

1807:   Level: advanced

1809: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1810: @*/
1811: PetscErrorCode TaoMonitorGradient(Tao tao, void *ctx)
1812: {
1813:   PetscViewer viewer = (PetscViewer)ctx;

1815:   PetscFunctionBegin;
1818:   PetscCall(VecView(tao->gradient, viewer));
1819:   PetscFunctionReturn(PETSC_SUCCESS);
1820: }

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

1825:   Collective

1827:   Input Parameters:
1828: + tao - the `Tao` context
1829: - ctx - `PetscViewer` context or `NULL`

1831:   Options Database Key:
1832: . -tao_monitor_step - view the step vector at each iteration

1834:   Level: advanced

1836: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1837: @*/
1838: PetscErrorCode TaoMonitorStep(Tao tao, void *ctx)
1839: {
1840:   PetscViewer viewer = (PetscViewer)ctx;

1842:   PetscFunctionBegin;
1845:   PetscCall(VecView(tao->stepdirection, viewer));
1846:   PetscFunctionReturn(PETSC_SUCCESS);
1847: }

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

1852:   Collective

1854:   Input Parameters:
1855: + tao - the `Tao` context
1856: - ctx - `TaoMonitorDraw` context

1858:   Options Database Key:
1859: . -tao_monitor_solution_draw - draw the solution at each iteration

1861:   Level: advanced

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

1867: .seealso: [](ch_tao), `Tao`, `TaoMonitorSolution()`, `TaoMonitorSet()`, `TaoMonitorGradientDraw()`, `TaoMonitorDrawCtxCreate()`,
1868:           `TaoMonitorDrawCtxDestroy()`
1869: @*/
1870: PetscErrorCode TaoMonitorSolutionDraw(Tao tao, void *ctx)
1871: {
1872:   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;

1874:   PetscFunctionBegin;
1876:   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(PETSC_SUCCESS);
1877:   PetscCall(VecView(tao->solution, ictx->viewer));
1878:   PetscFunctionReturn(PETSC_SUCCESS);
1879: }

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

1884:   Collective

1886:   Input Parameters:
1887: + tao - the `Tao` context
1888: - ctx - `PetscViewer` context

1890:   Options Database Key:
1891: . -tao_monitor_gradient_draw - draw the gradient at each iteration

1893:   Level: advanced

1895: .seealso: [](ch_tao), `Tao`, `TaoMonitorGradient()`, `TaoMonitorSet()`, `TaoMonitorSolutionDraw()`
1896: @*/
1897: PetscErrorCode TaoMonitorGradientDraw(Tao tao, void *ctx)
1898: {
1899:   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;

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

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

1911:   Collective

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

1917:   Options Database Key:
1918: . -tao_monitor_step_draw - draw the step direction at each iteration

1920:   Level: advanced

1922: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorSolutionDraw`
1923: @*/
1924: PetscErrorCode TaoMonitorStepDraw(Tao tao, void *ctx)
1925: {
1926:   PetscViewer viewer = (PetscViewer)ctx;

1928:   PetscFunctionBegin;
1931:   PetscCall(VecView(tao->stepdirection, viewer));
1932:   PetscFunctionReturn(PETSC_SUCCESS);
1933: }

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

1938:   Collective

1940:   Input Parameters:
1941: + tao - the `Tao` context
1942: - ctx - the `PetscViewer` context or `NULL`

1944:   Options Database Key:
1945: . -tao_monitor_ls_residual - view the residual at each iteration

1947:   Level: advanced

1949: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1950: @*/
1951: PetscErrorCode TaoMonitorResidual(Tao tao, void *ctx)
1952: {
1953:   PetscViewer viewer = (PetscViewer)ctx;

1955:   PetscFunctionBegin;
1958:   PetscCall(VecView(tao->ls_res, viewer));
1959:   PetscFunctionReturn(PETSC_SUCCESS);
1960: }

1962: /*@
1963:   TaoDefaultConvergenceTest - Determines whether the solver should continue iterating
1964:   or terminate.

1966:   Collective

1968:   Input Parameters:
1969: + tao   - the `Tao` context
1970: - dummy - unused dummy context

1972:   Level: developer

1974:   Notes:
1975:   This routine checks the residual in the optimality conditions, the
1976:   relative residual in the optimity conditions, the number of function
1977:   evaluations, and the function value to test convergence.  Some
1978:   solvers may use different convergence routines.

1980: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoGetConvergedReason()`, `TaoSetConvergedReason()`
1981: @*/
1982: PetscErrorCode TaoDefaultConvergenceTest(Tao tao, void *dummy)
1983: {
1984:   PetscInt           niter = tao->niter, nfuncs = PetscMax(tao->nfuncs, tao->nfuncgrads);
1985:   PetscInt           max_funcs = tao->max_funcs;
1986:   PetscReal          gnorm = tao->residual, gnorm0 = tao->gnorm0;
1987:   PetscReal          f = tao->fc, steptol = tao->steptol, trradius = tao->step;
1988:   PetscReal          gatol = tao->gatol, grtol = tao->grtol, gttol = tao->gttol;
1989:   PetscReal          catol = tao->catol, crtol = tao->crtol;
1990:   PetscReal          fmin = tao->fmin, cnorm = tao->cnorm;
1991:   TaoConvergedReason reason = tao->reason;

1993:   PetscFunctionBegin;
1995:   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);

1997:   if (PetscIsInfOrNanReal(f)) {
1998:     PetscCall(PetscInfo(tao, "Failed to converged, function value is Inf or NaN\n"));
1999:     reason = TAO_DIVERGED_NAN;
2000:   } else if (f <= fmin && cnorm <= catol) {
2001:     PetscCall(PetscInfo(tao, "Converged due to function value %g < minimum function value %g\n", (double)f, (double)fmin));
2002:     reason = TAO_CONVERGED_MINF;
2003:   } else if (gnorm <= gatol && cnorm <= catol) {
2004:     PetscCall(PetscInfo(tao, "Converged due to residual norm ||g(X)||=%g < %g\n", (double)gnorm, (double)gatol));
2005:     reason = TAO_CONVERGED_GATOL;
2006:   } else if (f != 0 && PetscAbsReal(gnorm / f) <= grtol && cnorm <= crtol) {
2007:     PetscCall(PetscInfo(tao, "Converged due to residual ||g(X)||/|f(X)| =%g < %g\n", (double)(gnorm / f), (double)grtol));
2008:     reason = TAO_CONVERGED_GRTOL;
2009:   } else if (gnorm0 != 0 && ((gttol == 0 && gnorm == 0) || gnorm / gnorm0 < gttol) && cnorm <= crtol) {
2010:     PetscCall(PetscInfo(tao, "Converged due to relative residual norm ||g(X)||/||g(X0)|| = %g < %g\n", (double)(gnorm / gnorm0), (double)gttol));
2011:     reason = TAO_CONVERGED_GTTOL;
2012:   } else if (max_funcs != PETSC_UNLIMITED && nfuncs > max_funcs) {
2013:     PetscCall(PetscInfo(tao, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", nfuncs, max_funcs));
2014:     reason = TAO_DIVERGED_MAXFCN;
2015:   } else if (tao->lsflag != 0) {
2016:     PetscCall(PetscInfo(tao, "Tao Line Search failure.\n"));
2017:     reason = TAO_DIVERGED_LS_FAILURE;
2018:   } else if (trradius < steptol && niter > 0) {
2019:     PetscCall(PetscInfo(tao, "Trust region/step size too small: %g < %g\n", (double)trradius, (double)steptol));
2020:     reason = TAO_CONVERGED_STEPTOL;
2021:   } else if (niter >= tao->max_it) {
2022:     PetscCall(PetscInfo(tao, "Exceeded maximum number of iterations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", niter, tao->max_it));
2023:     reason = TAO_DIVERGED_MAXITS;
2024:   } else {
2025:     reason = TAO_CONTINUE_ITERATING;
2026:   }
2027:   tao->reason = reason;
2028:   PetscFunctionReturn(PETSC_SUCCESS);
2029: }

2031: /*@
2032:   TaoSetOptionsPrefix - Sets the prefix used for searching for all
2033:   Tao options in the database.

2035:   Logically Collective

2037:   Input Parameters:
2038: + tao - the `Tao` context
2039: - p   - the prefix string to prepend to all Tao option requests

2041:   Level: advanced

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

2047:   For example, to distinguish between the runtime options for two
2048:   different Tao solvers, one could call
2049: .vb
2050:       TaoSetOptionsPrefix(tao1,"sys1_")
2051:       TaoSetOptionsPrefix(tao2,"sys2_")
2052: .ve

2054:   This would enable use of different options for each system, such as
2055: .vb
2056:       -sys1_tao_method blmvm -sys1_tao_grtol 1.e-3
2057:       -sys2_tao_method lmvm  -sys2_tao_grtol 1.e-4
2058: .ve

2060: .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoAppendOptionsPrefix()`, `TaoGetOptionsPrefix()`
2061: @*/
2062: PetscErrorCode TaoSetOptionsPrefix(Tao tao, const char p[])
2063: {
2064:   PetscFunctionBegin;
2066:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tao, p));
2067:   if (tao->linesearch) PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, p));
2068:   if (tao->ksp) PetscCall(KSPSetOptionsPrefix(tao->ksp, p));
2069:   PetscFunctionReturn(PETSC_SUCCESS);
2070: }

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

2075:   Logically Collective

2077:   Input Parameters:
2078: + tao - the `Tao` solver context
2079: - p   - the prefix string to prepend to all `Tao` option requests

2081:   Level: advanced

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

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

2099: /*@
2100:   TaoGetOptionsPrefix - Gets the prefix used for searching for all
2101:   Tao options in the database

2103:   Not Collective

2105:   Input Parameter:
2106: . tao - the `Tao` context

2108:   Output Parameter:
2109: . p - pointer to the prefix string used is returned

2111:   Fortran Notes:
2112:   Pass in a string 'prefix' of sufficient length to hold the prefix.

2114:   Level: advanced

2116: .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoSetOptionsPrefix()`, `TaoAppendOptionsPrefix()`
2117: @*/
2118: PetscErrorCode TaoGetOptionsPrefix(Tao tao, const char *p[])
2119: {
2120:   PetscFunctionBegin;
2122:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)tao, p));
2123:   PetscFunctionReturn(PETSC_SUCCESS);
2124: }

2126: /*@
2127:   TaoSetType - Sets the `TaoType` for the minimization solver.

2129:   Collective

2131:   Input Parameters:
2132: + tao  - the `Tao` solver context
2133: - type - a known method

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

2139:   Level: intermediate

2141: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoGetType()`, `TaoType`
2142: @*/
2143: PetscErrorCode TaoSetType(Tao tao, TaoType type)
2144: {
2145:   PetscErrorCode (*create_xxx)(Tao);
2146:   PetscBool issame;

2148:   PetscFunctionBegin;

2151:   PetscCall(PetscObjectTypeCompare((PetscObject)tao, type, &issame));
2152:   if (issame) PetscFunctionReturn(PETSC_SUCCESS);

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

2157:   /* Destroy the existing solver information */
2158:   PetscTryTypeMethod(tao, destroy);
2159:   PetscCall(KSPDestroy(&tao->ksp));
2160:   PetscCall(TaoLineSearchDestroy(&tao->linesearch));
2161:   tao->ops->setup          = NULL;
2162:   tao->ops->solve          = NULL;
2163:   tao->ops->view           = NULL;
2164:   tao->ops->setfromoptions = NULL;
2165:   tao->ops->destroy        = NULL;

2167:   tao->setupcalled = PETSC_FALSE;

2169:   PetscCall((*create_xxx)(tao));
2170:   PetscCall(PetscObjectChangeTypeName((PetscObject)tao, type));
2171:   PetscFunctionReturn(PETSC_SUCCESS);
2172: }

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

2177:   Not Collective, No Fortran Support

2179:   Input Parameters:
2180: + sname - name of a new user-defined solver
2181: - func  - routine to Create method context

2183:   Example Usage:
2184: .vb
2185:    TaoRegister("my_solver", MySolverCreate);
2186: .ve

2188:   Then, your solver can be chosen with the procedural interface via
2189: $     TaoSetType(tao, "my_solver")
2190:   or at runtime via the option
2191: $     -tao_type my_solver

2193:   Level: advanced

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

2198: .seealso: [](ch_tao), `Tao`, `TaoSetType()`, `TaoRegisterAll()`, `TaoRegisterDestroy()`
2199: @*/
2200: PetscErrorCode TaoRegister(const char sname[], PetscErrorCode (*func)(Tao))
2201: {
2202:   PetscFunctionBegin;
2203:   PetscCall(TaoInitializePackage());
2204:   PetscCall(PetscFunctionListAdd(&TaoList, sname, func));
2205:   PetscFunctionReturn(PETSC_SUCCESS);
2206: }

2208: /*@C
2209:   TaoRegisterDestroy - Frees the list of minimization solvers that were
2210:   registered by `TaoRegister()`.

2212:   Not Collective

2214:   Level: advanced

2216: .seealso: [](ch_tao), `Tao`, `TaoRegisterAll()`, `TaoRegister()`
2217: @*/
2218: PetscErrorCode TaoRegisterDestroy(void)
2219: {
2220:   PetscFunctionBegin;
2221:   PetscCall(PetscFunctionListDestroy(&TaoList));
2222:   TaoRegisterAllCalled = PETSC_FALSE;
2223:   PetscFunctionReturn(PETSC_SUCCESS);
2224: }

2226: /*@
2227:   TaoGetIterationNumber - Gets the number of `TaoSolve()` iterations completed
2228:   at this time.

2230:   Not Collective

2232:   Input Parameter:
2233: . tao - the `Tao` context

2235:   Output Parameter:
2236: . iter - iteration number

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

2241:   Level: intermediate

2243: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`, `TaoGetResidualNorm()`, `TaoGetObjective()`
2244: @*/
2245: PetscErrorCode TaoGetIterationNumber(Tao tao, PetscInt *iter)
2246: {
2247:   PetscFunctionBegin;
2249:   PetscAssertPointer(iter, 2);
2250:   *iter = tao->niter;
2251:   PetscFunctionReturn(PETSC_SUCCESS);
2252: }

2254: /*@
2255:   TaoGetResidualNorm - Gets the current value of the norm of the residual (gradient)
2256:   at this time.

2258:   Not Collective

2260:   Input Parameter:
2261: . tao - the `Tao` context

2263:   Output Parameter:
2264: . value - the current value

2266:   Level: intermediate

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

2272: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`, `TaoGetIterationNumber()`, `TaoGetObjective()`
2273: @*/
2274: PetscErrorCode TaoGetResidualNorm(Tao tao, PetscReal *value)
2275: {
2276:   PetscFunctionBegin;
2278:   PetscAssertPointer(value, 2);
2279:   *value = tao->residual;
2280:   PetscFunctionReturn(PETSC_SUCCESS);
2281: }

2283: /*@
2284:   TaoSetIterationNumber - Sets the current iteration number.

2286:   Logically Collective

2288:   Input Parameters:
2289: + tao  - the `Tao` context
2290: - iter - iteration number

2292:   Level: developer

2294: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
2295: @*/
2296: PetscErrorCode TaoSetIterationNumber(Tao tao, PetscInt iter)
2297: {
2298:   PetscFunctionBegin;
2301:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2302:   tao->niter = iter;
2303:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2304:   PetscFunctionReturn(PETSC_SUCCESS);
2305: }

2307: /*@
2308:   TaoGetTotalIterationNumber - Gets the total number of `TaoSolve()` iterations
2309:   completed. This number keeps accumulating if multiple solves
2310:   are called with the `Tao` object.

2312:   Not Collective

2314:   Input Parameter:
2315: . tao - the `Tao` context

2317:   Output Parameter:
2318: . iter - number of iterations

2320:   Level: intermediate

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

2326: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
2327: @*/
2328: PetscErrorCode TaoGetTotalIterationNumber(Tao tao, PetscInt *iter)
2329: {
2330:   PetscFunctionBegin;
2332:   PetscAssertPointer(iter, 2);
2333:   *iter = tao->ntotalits;
2334:   PetscFunctionReturn(PETSC_SUCCESS);
2335: }

2337: /*@
2338:   TaoSetTotalIterationNumber - Sets the current total iteration number.

2340:   Logically Collective

2342:   Input Parameters:
2343: + tao  - the `Tao` context
2344: - iter - the iteration number

2346:   Level: developer

2348: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
2349: @*/
2350: PetscErrorCode TaoSetTotalIterationNumber(Tao tao, PetscInt iter)
2351: {
2352:   PetscFunctionBegin;
2355:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2356:   tao->ntotalits = iter;
2357:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2358:   PetscFunctionReturn(PETSC_SUCCESS);
2359: }

2361: /*@
2362:   TaoSetConvergedReason - Sets the termination flag on a `Tao` object

2364:   Logically Collective

2366:   Input Parameters:
2367: + tao    - the `Tao` context
2368: - reason - the `TaoConvergedReason`

2370:   Level: intermediate

2372: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`
2373: @*/
2374: PetscErrorCode TaoSetConvergedReason(Tao tao, TaoConvergedReason reason)
2375: {
2376:   PetscFunctionBegin;
2379:   tao->reason = reason;
2380:   PetscFunctionReturn(PETSC_SUCCESS);
2381: }

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

2386:   Not Collective

2388:   Input Parameter:
2389: . tao - the `Tao` solver context

2391:   Output Parameter:
2392: . reason - value of `TaoConvergedReason`

2394:   Level: intermediate

2396: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetConvergenceTest()`, `TaoSetTolerances()`
2397: @*/
2398: PetscErrorCode TaoGetConvergedReason(Tao tao, TaoConvergedReason *reason)
2399: {
2400:   PetscFunctionBegin;
2402:   PetscAssertPointer(reason, 2);
2403:   *reason = tao->reason;
2404:   PetscFunctionReturn(PETSC_SUCCESS);
2405: }

2407: /*@
2408:   TaoGetSolutionStatus - Get the current iterate, objective value,
2409:   residual, infeasibility, and termination from a `Tao` object

2411:   Not Collective

2413:   Input Parameter:
2414: . tao - the `Tao` context

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

2424:   Level: intermediate

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

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

2431: .seealso: [](ch_tao), `TaoMonitor()`, `TaoGetConvergedReason()`
2432: @*/
2433: PetscErrorCode TaoGetSolutionStatus(Tao tao, PetscInt *its, PetscReal *f, PetscReal *gnorm, PetscReal *cnorm, PetscReal *xdiff, TaoConvergedReason *reason)
2434: {
2435:   PetscFunctionBegin;
2437:   if (its) *its = tao->niter;
2438:   if (f) *f = tao->fc;
2439:   if (gnorm) *gnorm = tao->residual;
2440:   if (cnorm) *cnorm = tao->cnorm;
2441:   if (reason) *reason = tao->reason;
2442:   if (xdiff) *xdiff = tao->step;
2443:   PetscFunctionReturn(PETSC_SUCCESS);
2444: }

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

2449:   Not Collective

2451:   Input Parameter:
2452: . tao - the `Tao` solver context

2454:   Output Parameter:
2455: . type - the `TaoType`

2457:   Level: intermediate

2459: .seealso: [](ch_tao), `Tao`, `TaoType`, `TaoSetType()`
2460: @*/
2461: PetscErrorCode TaoGetType(Tao tao, TaoType *type)
2462: {
2463:   PetscFunctionBegin;
2465:   PetscAssertPointer(type, 2);
2466:   *type = ((PetscObject)tao)->type_name;
2467:   PetscFunctionReturn(PETSC_SUCCESS);
2468: }

2470: /*@C
2471:   TaoMonitor - Monitor the solver and the current solution.  This
2472:   routine will record the iteration number and residual statistics,
2473:   and call any monitors specified by the user.

2475:   Input Parameters:
2476: + tao        - the `Tao` context
2477: . its        - the current iterate number (>=0)
2478: . f          - the current objective function value
2479: . res        - the gradient norm, square root of the duality gap, or other measure indicating distance from optimality.  This measure will be recorded and
2480:           used for some termination tests.
2481: . cnorm      - the infeasibility of the current solution with regard to the constraints.
2482: - steplength - multiple of the step direction added to the previous iterate.

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

2487:   Level: developer

2489: .seealso: [](ch_tao), `Tao`, `TaoGetConvergedReason()`, `TaoMonitorDefault()`, `TaoMonitorSet()`
2490: @*/
2491: PetscErrorCode TaoMonitor(Tao tao, PetscInt its, PetscReal f, PetscReal res, PetscReal cnorm, PetscReal steplength)
2492: {
2493:   PetscInt i;

2495:   PetscFunctionBegin;
2497:   tao->fc       = f;
2498:   tao->residual = res;
2499:   tao->cnorm    = cnorm;
2500:   tao->step     = steplength;
2501:   if (!its) {
2502:     tao->cnorm0 = cnorm;
2503:     tao->gnorm0 = res;
2504:   }
2505:   PetscCall(VecLockReadPush(tao->solution));
2506:   for (i = 0; i < tao->numbermonitors; i++) PetscCall((*tao->monitor[i])(tao, tao->monitorcontext[i]));
2507:   PetscCall(VecLockReadPop(tao->solution));
2508:   PetscFunctionReturn(PETSC_SUCCESS);
2509: }

2511: /*@
2512:   TaoSetConvergenceHistory - Sets the array used to hold the convergence history.

2514:   Logically Collective

2516:   Input Parameters:
2517: + tao   - the `Tao` solver context
2518: . obj   - array to hold objective value history
2519: . resid - array to hold residual history
2520: . cnorm - array to hold constraint violation history
2521: . lits  - integer array holds the number of linear iterations for each Tao iteration
2522: . na    - size of `obj`, `resid`, and `cnorm`
2523: - reset - `PETSC_TRUE` indicates each new minimization resets the history counter to zero,
2524:            else it continues storing new values for new minimizations after the old ones

2526:   Level: intermediate

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

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

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

2541: .seealso: [](ch_tao), `TaoGetConvergenceHistory()`
2542: @*/
2543: PetscErrorCode TaoSetConvergenceHistory(Tao tao, PetscReal obj[], PetscReal resid[], PetscReal cnorm[], PetscInt lits[], PetscInt na, PetscBool reset)
2544: {
2545:   PetscFunctionBegin;
2547:   if (obj) PetscAssertPointer(obj, 2);
2548:   if (resid) PetscAssertPointer(resid, 3);
2549:   if (cnorm) PetscAssertPointer(cnorm, 4);
2550:   if (lits) PetscAssertPointer(lits, 5);

2552:   if (na == PETSC_DECIDE || na == PETSC_CURRENT) na = 1000;
2553:   if (!obj && !resid && !cnorm && !lits) {
2554:     PetscCall(PetscCalloc4(na, &obj, na, &resid, na, &cnorm, na, &lits));
2555:     tao->hist_malloc = PETSC_TRUE;
2556:   }

2558:   tao->hist_obj   = obj;
2559:   tao->hist_resid = resid;
2560:   tao->hist_cnorm = cnorm;
2561:   tao->hist_lits  = lits;
2562:   tao->hist_max   = na;
2563:   tao->hist_reset = reset;
2564:   tao->hist_len   = 0;
2565:   PetscFunctionReturn(PETSC_SUCCESS);
2566: }

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

2571:   Collective

2573:   Input Parameter:
2574: . tao - the `Tao` context

2576:   Output Parameters:
2577: + obj   - array used to hold objective value history
2578: . resid - array used to hold residual history
2579: . cnorm - array used to hold constraint violation history
2580: . lits  - integer array used to hold linear solver iteration count
2581: - nhist - size of `obj`, `resid`, `cnorm`, and `lits`

2583:   Level: advanced

2585:   Notes:
2586:   This routine must be preceded by calls to `TaoSetConvergenceHistory()`
2587:   and `TaoSolve()`, otherwise it returns useless information.

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

2593:   Fortran Notes:
2594:   The calling sequence is
2595: .vb
2596:    call TaoGetConvergenceHistory(Tao tao, PetscInt nhist, PetscErrorCode ierr)
2597: .ve
2598:   In other words this gets the current number of entries in the history. Access the history through the array you passed to `TaoSetConvergenceHistory()`

2600: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetConvergenceHistory()`
2601: @*/
2602: PetscErrorCode TaoGetConvergenceHistory(Tao tao, PetscReal **obj, PetscReal **resid, PetscReal **cnorm, PetscInt **lits, PetscInt *nhist)
2603: {
2604:   PetscFunctionBegin;
2606:   if (obj) *obj = tao->hist_obj;
2607:   if (cnorm) *cnorm = tao->hist_cnorm;
2608:   if (resid) *resid = tao->hist_resid;
2609:   if (lits) *lits = tao->hist_lits;
2610:   if (nhist) *nhist = tao->hist_len;
2611:   PetscFunctionReturn(PETSC_SUCCESS);
2612: }

2614: /*@
2615:   TaoSetApplicationContext - Sets the optional user-defined context for a `Tao` solver.

2617:   Logically Collective

2619:   Input Parameters:
2620: + tao - the `Tao` context
2621: - ctx - optional user context

2623:   Level: intermediate

2625: .seealso: [](ch_tao), `Tao`, `TaoGetApplicationContext()`
2626: @*/
2627: PetscErrorCode TaoSetApplicationContext(Tao tao, void *ctx)
2628: {
2629:   PetscFunctionBegin;
2631:   tao->ctx = ctx;
2632:   PetscFunctionReturn(PETSC_SUCCESS);
2633: }

2635: /*@
2636:   TaoGetApplicationContext - Gets the user-defined context for a `Tao` solver

2638:   Not Collective

2640:   Input Parameter:
2641: . tao - the `Tao` context

2643:   Output Parameter:
2644: . ctx - user context

2646:   Level: intermediate

2648: .seealso: [](ch_tao), `Tao`, `TaoSetApplicationContext()`
2649: @*/
2650: PetscErrorCode TaoGetApplicationContext(Tao tao, void *ctx)
2651: {
2652:   PetscFunctionBegin;
2654:   PetscAssertPointer(ctx, 2);
2655:   *(void **)ctx = tao->ctx;
2656:   PetscFunctionReturn(PETSC_SUCCESS);
2657: }

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

2662:   Collective

2664:   Input Parameters:
2665: + tao - the `Tao` context
2666: - M   - matrix that defines the norm

2668:   Level: beginner

2670: .seealso: [](ch_tao), `Tao`, `TaoGetGradientNorm()`, `TaoGradientNorm()`
2671: @*/
2672: PetscErrorCode TaoSetGradientNorm(Tao tao, Mat M)
2673: {
2674:   PetscFunctionBegin;
2677:   PetscCall(PetscObjectReference((PetscObject)M));
2678:   PetscCall(MatDestroy(&tao->gradient_norm));
2679:   PetscCall(VecDestroy(&tao->gradient_norm_tmp));
2680:   tao->gradient_norm = M;
2681:   PetscCall(MatCreateVecs(M, NULL, &tao->gradient_norm_tmp));
2682:   PetscFunctionReturn(PETSC_SUCCESS);
2683: }

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

2688:   Not Collective

2690:   Input Parameter:
2691: . tao - the `Tao` context

2693:   Output Parameter:
2694: . M - gradient norm

2696:   Level: beginner

2698: .seealso: [](ch_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGradientNorm()`
2699: @*/
2700: PetscErrorCode TaoGetGradientNorm(Tao tao, Mat *M)
2701: {
2702:   PetscFunctionBegin;
2704:   PetscAssertPointer(M, 2);
2705:   *M = tao->gradient_norm;
2706:   PetscFunctionReturn(PETSC_SUCCESS);
2707: }

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

2712:   Collective

2714:   Input Parameters:
2715: + tao      - the `Tao` context
2716: . gradient - the gradient
2717: - type     - the norm type

2719:   Output Parameter:
2720: . gnorm - the gradient norm

2722:   Level: advanced

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

2727:   Developer Notes:
2728:   Should be named `TaoComputeGradientNorm()`.

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

2733: .seealso: [](ch_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGetGradientNorm()`
2734: @*/
2735: PetscErrorCode TaoGradientNorm(Tao tao, Vec gradient, NormType type, PetscReal *gnorm)
2736: {
2737:   PetscFunctionBegin;
2741:   PetscAssertPointer(gnorm, 4);
2742:   if (tao->gradient_norm) {
2743:     PetscScalar gnorms;

2745:     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.");
2746:     PetscCall(MatMult(tao->gradient_norm, gradient, tao->gradient_norm_tmp));
2747:     PetscCall(VecDot(gradient, tao->gradient_norm_tmp, &gnorms));
2748:     *gnorm = PetscRealPart(PetscSqrtScalar(gnorms));
2749:   } else {
2750:     PetscCall(VecNorm(gradient, type, gnorm));
2751:   }
2752:   PetscFunctionReturn(PETSC_SUCCESS);
2753: }

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

2758:   Collective

2760:   Input Parameters:
2761: + comm     - the communicator to share the context
2762: . host     - the name of the X Windows host that will display the monitor
2763: . label    - the label to put at the top of the display window
2764: . x        - the horizontal coordinate of the lower left corner of the window to open
2765: . y        - the vertical coordinate of the lower left corner of the window to open
2766: . m        - the width of the window
2767: . n        - the height of the window
2768: - howoften - how many `Tao` iterations between displaying the monitor information

2770:   Output Parameter:
2771: . ctx - the monitor context

2773:   Options Database Keys:
2774: + -tao_monitor_solution_draw - use `TaoMonitorSolutionDraw()` to monitor the solution
2775: - -tao_draw_solution_initial - show initial guess as well as current solution

2777:   Level: intermediate

2779:   Note:
2780:   The context this creates, along with `TaoMonitorSolutionDraw()`, and `TaoMonitorDrawCtxDestroy()`
2781:   are passed to `TaoMonitorSet()`.

2783: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorDrawCtx()`
2784: @*/
2785: PetscErrorCode TaoMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TaoMonitorDrawCtx *ctx)
2786: {
2787:   PetscFunctionBegin;
2788:   PetscCall(PetscNew(ctx));
2789:   PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
2790:   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
2791:   (*ctx)->howoften = howoften;
2792:   PetscFunctionReturn(PETSC_SUCCESS);
2793: }

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

2798:   Collective

2800:   Input Parameter:
2801: . ictx - the monitor context

2803:   Level: intermediate

2805:   Note:
2806:   This is passed to `TaoMonitorSet()` as the final argument, along with `TaoMonitorSolutionDraw()`, and the context
2807:   obtained with `TaoMonitorDrawCtxCreate()`.

2809: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorSolutionDraw()`
2810: @*/
2811: PetscErrorCode TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx *ictx)
2812: {
2813:   PetscFunctionBegin;
2814:   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
2815:   PetscCall(PetscFree(*ictx));
2816:   PetscFunctionReturn(PETSC_SUCCESS);
2817: }