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:   TaoCreate - Creates a Tao solver

 75:   Collective

 77:   Input Parameter:
 78: . comm - MPI communicator

 80:   Output Parameter:
 81: . newtao - the new `Tao` context

 83:   Options Database Key:
 84: . -tao_type - select which method Tao should use

 86:   Level: beginner

 88: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoDestroy()`, `TaoSetFromOptions()`, `TaoSetType()`
 89: @*/
 90: PetscErrorCode TaoCreate(MPI_Comm comm, Tao *newtao)
 91: {
 92:   Tao tao;

 94:   PetscFunctionBegin;
 95:   PetscAssertPointer(newtao, 2);
 96:   PetscCall(TaoInitializePackage());
 97:   PetscCall(TaoLineSearchInitializePackage());
 98:   PetscCall(PetscHeaderCreate(tao, TAO_CLASSID, "Tao", "Optimization solver", "Tao", comm, TaoDestroy, TaoView));

100:   /* Set non-NULL defaults */
101:   tao->ops->convergencetest = TaoDefaultConvergenceTest;

103:   tao->max_it    = 10000;
104:   tao->max_funcs = -1;
105: #if defined(PETSC_USE_REAL_SINGLE)
106:   tao->gatol = 1e-5;
107:   tao->grtol = 1e-5;
108:   tao->crtol = 1e-5;
109:   tao->catol = 1e-5;
110: #else
111:   tao->gatol = 1e-8;
112:   tao->grtol = 1e-8;
113:   tao->crtol = 1e-8;
114:   tao->catol = 1e-8;
115: #endif
116:   tao->gttol   = 0.0;
117:   tao->steptol = 0.0;
118:   tao->trust0  = PETSC_INFINITY;
119:   tao->fmin    = PETSC_NINFINITY;

121:   tao->hist_reset = PETSC_TRUE;

123:   PetscCall(TaoResetStatistics(tao));
124:   *newtao = tao;
125:   PetscFunctionReturn(PETSC_SUCCESS);
126: }

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

131:   Collective

133:   Input Parameter:
134: . tao - the `Tao` context

136:   Level: beginner

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

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

144: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSetObjective()`, `TaoSetGradient()`, `TaoSetHessian()`, `TaoGetConvergedReason()`, `TaoSetUp()`
145:  @*/
146: PetscErrorCode TaoSolve(Tao tao)
147: {
148:   static PetscBool set = PETSC_FALSE;

150:   PetscFunctionBegin;
152:   PetscCall(PetscCitationsRegister("@TechReport{tao-user-ref,\n"
153:                                    "title   = {Toolkit for Advanced Optimization (TAO) Users Manual},\n"
154:                                    "author  = {Todd Munson and Jason Sarich and Stefan Wild and Steve Benson and Lois Curfman McInnes},\n"
155:                                    "Institution = {Argonne National Laboratory},\n"
156:                                    "Year   = 2014,\n"
157:                                    "Number = {ANL/MCS-TM-322 - Revision 3.5},\n"
158:                                    "url    = {https://www.mcs.anl.gov/research/projects/tao/}\n}\n",
159:                                    &set));
160:   tao->header_printed = PETSC_FALSE;
161:   PetscCall(TaoSetUp(tao));
162:   PetscCall(TaoResetStatistics(tao));
163:   if (tao->linesearch) PetscCall(TaoLineSearchReset(tao->linesearch));

165:   PetscCall(PetscLogEventBegin(TAO_Solve, tao, 0, 0, 0));
166:   PetscTryTypeMethod(tao, solve);
167:   PetscCall(PetscLogEventEnd(TAO_Solve, tao, 0, 0, 0));

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

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

173:   if (tao->printreason) {
174:     PetscViewer viewer = PETSC_VIEWER_STDOUT_(((PetscObject)tao)->comm);
175:     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)tao)->tablevel));
176:     if (tao->reason > 0) {
177:       PetscCall(PetscViewerASCIIPrintf(viewer, "  TAO %s solve converged due to %s iterations %" PetscInt_FMT "\n", ((PetscObject)tao)->prefix ? ((PetscObject)tao)->prefix : "", TaoConvergedReasons[tao->reason], tao->niter));
178:     } else {
179:       PetscCall(PetscViewerASCIIPrintf(viewer, "  TAO %s solve did not converge due to %s iteration %" PetscInt_FMT "\n", ((PetscObject)tao)->prefix ? ((PetscObject)tao)->prefix : "", TaoConvergedReasons[tao->reason], tao->niter));
180:     }
181:     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)tao)->tablevel));
182:   }
183:   PetscCall(TaoViewFromOptions(tao, NULL, "-tao_view"));
184:   PetscFunctionReturn(PETSC_SUCCESS);
185: }

187: /*@
188:   TaoSetUp - Sets up the internal data structures for the later use
189:   of a Tao solver

191:   Collective

193:   Input Parameter:
194: . tao - the `Tao` context

196:   Level: advanced

198:   Note:
199:   The user will not need to explicitly call `TaoSetUp()`, as it will
200:   automatically be called in `TaoSolve()`.  However, if the user
201:   desires to call it explicitly, it should come after `TaoCreate()`
202:   and any TaoSetSomething() routines, but before `TaoSolve()`.

204: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
205: @*/
206: PetscErrorCode TaoSetUp(Tao tao)
207: {
208:   PetscFunctionBegin;
210:   if (tao->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
211:   PetscCall(TaoSetUpEW_Private(tao));
212:   PetscCheck(tao->solution, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoSetSolution");
213:   PetscTryTypeMethod(tao, setup);
214:   tao->setupcalled = PETSC_TRUE;
215:   PetscFunctionReturn(PETSC_SUCCESS);
216: }

218: /*@C
219:   TaoDestroy - Destroys the `Tao` context that was created with `TaoCreate()`

221:   Collective

223:   Input Parameter:
224: . tao - the `Tao` context

226:   Level: beginner

228: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
229: @*/
230: PetscErrorCode TaoDestroy(Tao *tao)
231: {
232:   PetscFunctionBegin;
233:   if (!*tao) PetscFunctionReturn(PETSC_SUCCESS);
235:   if (--((PetscObject)*tao)->refct > 0) {
236:     *tao = NULL;
237:     PetscFunctionReturn(PETSC_SUCCESS);
238:   }

240:   PetscTryTypeMethod(*tao, destroy);
241:   PetscCall(KSPDestroy(&(*tao)->ksp));
242:   PetscCall(SNESDestroy(&(*tao)->snes_ewdummy));
243:   PetscCall(TaoLineSearchDestroy(&(*tao)->linesearch));

245:   if ((*tao)->ops->convergencedestroy) {
246:     PetscCall((*(*tao)->ops->convergencedestroy)((*tao)->cnvP));
247:     if ((*tao)->jacobian_state_inv) PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
248:   }
249:   PetscCall(VecDestroy(&(*tao)->solution));
250:   PetscCall(VecDestroy(&(*tao)->gradient));
251:   PetscCall(VecDestroy(&(*tao)->ls_res));

253:   if ((*tao)->gradient_norm) {
254:     PetscCall(PetscObjectDereference((PetscObject)(*tao)->gradient_norm));
255:     PetscCall(VecDestroy(&(*tao)->gradient_norm_tmp));
256:   }

258:   PetscCall(VecDestroy(&(*tao)->XL));
259:   PetscCall(VecDestroy(&(*tao)->XU));
260:   PetscCall(VecDestroy(&(*tao)->IL));
261:   PetscCall(VecDestroy(&(*tao)->IU));
262:   PetscCall(VecDestroy(&(*tao)->DE));
263:   PetscCall(VecDestroy(&(*tao)->DI));
264:   PetscCall(VecDestroy(&(*tao)->constraints));
265:   PetscCall(VecDestroy(&(*tao)->constraints_equality));
266:   PetscCall(VecDestroy(&(*tao)->constraints_inequality));
267:   PetscCall(VecDestroy(&(*tao)->stepdirection));
268:   PetscCall(MatDestroy(&(*tao)->hessian_pre));
269:   PetscCall(MatDestroy(&(*tao)->hessian));
270:   PetscCall(MatDestroy(&(*tao)->ls_jac));
271:   PetscCall(MatDestroy(&(*tao)->ls_jac_pre));
272:   PetscCall(MatDestroy(&(*tao)->jacobian_pre));
273:   PetscCall(MatDestroy(&(*tao)->jacobian));
274:   PetscCall(MatDestroy(&(*tao)->jacobian_state_pre));
275:   PetscCall(MatDestroy(&(*tao)->jacobian_state));
276:   PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
277:   PetscCall(MatDestroy(&(*tao)->jacobian_design));
278:   PetscCall(MatDestroy(&(*tao)->jacobian_equality));
279:   PetscCall(MatDestroy(&(*tao)->jacobian_equality_pre));
280:   PetscCall(MatDestroy(&(*tao)->jacobian_inequality));
281:   PetscCall(MatDestroy(&(*tao)->jacobian_inequality_pre));
282:   PetscCall(ISDestroy(&(*tao)->state_is));
283:   PetscCall(ISDestroy(&(*tao)->design_is));
284:   PetscCall(VecDestroy(&(*tao)->res_weights_v));
285:   PetscCall(TaoMonitorCancel(*tao));
286:   if ((*tao)->hist_malloc) PetscCall(PetscFree4((*tao)->hist_obj, (*tao)->hist_resid, (*tao)->hist_cnorm, (*tao)->hist_lits));
287:   if ((*tao)->res_weights_n) {
288:     PetscCall(PetscFree((*tao)->res_weights_rows));
289:     PetscCall(PetscFree((*tao)->res_weights_cols));
290:     PetscCall(PetscFree((*tao)->res_weights_w));
291:   }
292:   PetscCall(PetscHeaderDestroy(tao));
293:   PetscFunctionReturn(PETSC_SUCCESS);
294: }

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

299:   Logically Collective

301:   Input Parameters:
302: + tao  - Tao context
303: - flag - `PETSC_TRUE` or `PETSC_FALSE`

305:   Level: advanced

307:   Note:
308:   See `SNESKSPSetUseEW()` for customization details.

310: .seealso: [](ch_tao), `Tao`, `SNESKSPSetUseEW()`
311: @*/
312: PetscErrorCode TaoKSPSetUseEW(Tao tao, PetscBool flag)
313: {
314:   PetscFunctionBegin;
317:   tao->ksp_ewconv = flag;
318:   PetscFunctionReturn(PETSC_SUCCESS);
319: }

321: /*@
322:   TaoSetFromOptions - Sets various Tao parameters from the options database

324:   Collective

326:   Input Parameter:
327: . tao - the `Tao` solver context

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

358:   Level: beginner

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

364: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
365: @*/
366: PetscErrorCode TaoSetFromOptions(Tao tao)
367: {
368:   TaoType     default_type = TAOLMVM;
369:   char        type[256], monfilename[PETSC_MAX_PATH_LEN];
370:   PetscViewer monviewer;
371:   PetscBool   flg, found;
372:   MPI_Comm    comm;

374:   PetscFunctionBegin;
376:   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));

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

380:   PetscObjectOptionsBegin((PetscObject)tao);
381:   /* Check for type from options */
382:   PetscCall(PetscOptionsFList("-tao_type", "Tao Solver type", "TaoSetType", TaoList, default_type, type, 256, &flg));
383:   if (flg) {
384:     PetscCall(TaoSetType(tao, type));
385:   } else if (!((PetscObject)tao)->type_name) {
386:     PetscCall(TaoSetType(tao, default_type));
387:   }

389:   /* Tao solvers do not set the prefix, set it here if not yet done
390:      We do it after SetType since solver may have been changed */
391:   if (tao->linesearch) {
392:     const char *prefix;
393:     PetscCall(TaoLineSearchGetOptionsPrefix(tao->linesearch, &prefix));
394:     if (!prefix) PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, ((PetscObject)tao)->prefix));
395:   }

397:   PetscCall(PetscOptionsReal("-tao_catol", "Stop if constraints violations within", "TaoSetConstraintTolerances", tao->catol, &tao->catol, &flg));
398:   if (flg) tao->catol_changed = PETSC_TRUE;
399:   PetscCall(PetscOptionsReal("-tao_crtol", "Stop if relative constraint violations within", "TaoSetConstraintTolerances", tao->crtol, &tao->crtol, &flg));
400:   if (flg) tao->crtol_changed = PETSC_TRUE;
401:   PetscCall(PetscOptionsReal("-tao_gatol", "Stop if norm of gradient less than", "TaoSetTolerances", tao->gatol, &tao->gatol, &flg));
402:   if (flg) tao->gatol_changed = PETSC_TRUE;
403:   PetscCall(PetscOptionsReal("-tao_grtol", "Stop if norm of gradient divided by the function value is less than", "TaoSetTolerances", tao->grtol, &tao->grtol, &flg));
404:   if (flg) tao->grtol_changed = PETSC_TRUE;
405:   PetscCall(PetscOptionsReal("-tao_gttol", "Stop if the norm of the gradient is less than the norm of the initial gradient times tol", "TaoSetTolerances", tao->gttol, &tao->gttol, &flg));
406:   if (flg) tao->gttol_changed = PETSC_TRUE;
407:   PetscCall(PetscOptionsInt("-tao_max_it", "Stop if iteration number exceeds", "TaoSetMaximumIterations", tao->max_it, &tao->max_it, &flg));
408:   if (flg) tao->max_it_changed = PETSC_TRUE;
409:   PetscCall(PetscOptionsInt("-tao_max_funcs", "Stop if number of function evaluations exceeds", "TaoSetMaximumFunctionEvaluations", tao->max_funcs, &tao->max_funcs, &flg));
410:   if (flg) tao->max_funcs_changed = PETSC_TRUE;
411:   PetscCall(PetscOptionsReal("-tao_fmin", "Stop if function less than", "TaoSetFunctionLowerBound", tao->fmin, &tao->fmin, &flg));
412:   if (flg) tao->fmin_changed = PETSC_TRUE;
413:   PetscCall(PetscOptionsReal("-tao_steptol", "Stop if step size or trust region radius less than", "", tao->steptol, &tao->steptol, &flg));
414:   if (flg) tao->steptol_changed = PETSC_TRUE;
415:   PetscCall(PetscOptionsReal("-tao_trust0", "Initial trust region radius", "TaoSetTrustRegionRadius", tao->trust0, &tao->trust0, &flg));
416:   if (flg) tao->trust0_changed = PETSC_TRUE;

418:   PetscCall(PetscOptionsDeprecated("-tao_solution_monitor", "-tao_monitor_solution", "3.21", NULL));
419:   PetscCall(PetscOptionsDeprecated("-tao_gradient_monitor", "-tao_monitor_gradient", "3.21", NULL));
420:   PetscCall(PetscOptionsDeprecated("-tao_stepdirection_monitor", "-tao_monitor_step", "3.21", NULL));
421:   PetscCall(PetscOptionsDeprecated("-tao_residual_monitor", "-tao_monitor_residual", "3.21", NULL));
422:   PetscCall(PetscOptionsDeprecated("-tao_smonitor", "-tao_monitor_short", "3.21", NULL));
423:   PetscCall(PetscOptionsDeprecated("-tao_cmonitor", "-tao_monitor_constraint_norm", "3.21", NULL));
424:   PetscCall(PetscOptionsDeprecated("-tao_gmonitor", "-tao_monitor_globalization", "3.21", NULL));
425:   PetscCall(PetscOptionsDeprecated("-tao_draw_solution", "-tao_monitor_solution_draw", "3.21", NULL));
426:   PetscCall(PetscOptionsDeprecated("-tao_draw_gradient", "-tao_monitor_gradient_draw", "3.21", NULL));
427:   PetscCall(PetscOptionsDeprecated("-tao_draw_step", "-tao_monitor_step_draw", "3.21", NULL));

429:   PetscCall(PetscOptionsString("-tao_monitor_solution", "View solution vector after each iteration", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
430:   if (flg) {
431:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
432:     PetscCall(TaoMonitorSet(tao, TaoMonitorSolution, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
433:   }

435:   PetscCall(PetscOptionsBool("-tao_converged_reason", "Print reason for Tao converged", "TaoSolve", tao->printreason, &tao->printreason, NULL));
436:   PetscCall(PetscOptionsString("-tao_monitor_gradient", "View gradient vector for each iteration", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
437:   if (flg) {
438:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
439:     PetscCall(TaoMonitorSet(tao, TaoMonitorGradient, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
440:   }

442:   PetscCall(PetscOptionsString("-tao_monitor_step", "View step vector after each iteration", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
443:   if (flg) {
444:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
445:     PetscCall(TaoMonitorSet(tao, TaoMonitorStep, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
446:   }

448:   PetscCall(PetscOptionsString("-tao_monitor_residual", "View least-squares residual vector after each iteration", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
449:   if (flg) {
450:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
451:     PetscCall(TaoMonitorSet(tao, TaoMonitorResidual, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
452:   }

454:   PetscCall(PetscOptionsString("-tao_monitor", "Use the default convergence monitor", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
455:   if (flg) {
456:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
457:     PetscCall(TaoMonitorSet(tao, TaoMonitorDefault, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
458:   }

460:   PetscCall(PetscOptionsString("-tao_monitor_globalization", "Use the convergence monitor with extra globalization info", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
461:   if (flg) {
462:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
463:     PetscCall(TaoMonitorSet(tao, TaoMonitorGlobalization, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
464:   }

466:   PetscCall(PetscOptionsString("-tao_monitor_short", "Use the short convergence monitor", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
467:   if (flg) {
468:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
469:     PetscCall(TaoMonitorSet(tao, TaoMonitorDefaultShort, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
470:   }

472:   PetscCall(PetscOptionsString("-tao_monitor_constraint_norm", "Use the default convergence monitor with constraint norm", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
473:   if (flg) {
474:     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
475:     PetscCall(TaoMonitorSet(tao, TaoMonitorConstraintNorm, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
476:   }

478:   flg = PETSC_FALSE;
479:   PetscCall(PetscOptionsDeprecated("-tao_cancelmonitors", "-tao_monitor_cancel", "3.21", NULL));
480:   PetscCall(PetscOptionsBool("-tao_monitor_cancel", "cancel all monitors and call any registered destroy routines", "TaoMonitorCancel", flg, &flg, NULL));
481:   if (flg) PetscCall(TaoMonitorCancel(tao));

483:   flg = PETSC_FALSE;
484:   PetscCall(PetscOptionsBool("-tao_monitor_solution_draw", "Plot solution vector at each iteration", "TaoMonitorSet", flg, &flg, NULL));
485:   if (flg) {
486:     TaoMonitorDrawCtx drawctx;
487:     PetscInt          howoften = 1;
488:     PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &drawctx));
489:     PetscCall(TaoMonitorSet(tao, TaoMonitorSolutionDraw, drawctx, (PetscErrorCode(*)(void **))TaoMonitorDrawCtxDestroy));
490:   }

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

496:   flg = PETSC_FALSE;
497:   PetscCall(PetscOptionsBool("-tao_monitor_gradient_draw", "plots gradient at each iteration", "TaoMonitorSet", flg, &flg, NULL));
498:   if (flg) {
499:     TaoMonitorDrawCtx drawctx;
500:     PetscInt          howoften = 1;
501:     PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &drawctx));
502:     PetscCall(TaoMonitorSet(tao, TaoMonitorGradientDraw, drawctx, (PetscErrorCode(*)(void **))TaoMonitorDrawCtxDestroy));
503:   }
504:   flg = PETSC_FALSE;
505:   PetscCall(PetscOptionsBool("-tao_fd_gradient", "compute gradient using finite differences", "TaoDefaultComputeGradient", flg, &flg, NULL));
506:   if (flg) PetscCall(TaoSetGradient(tao, NULL, TaoDefaultComputeGradient, NULL));
507:   flg = PETSC_FALSE;
508:   PetscCall(PetscOptionsBool("-tao_fd_hessian", "compute Hessian using finite differences", "TaoDefaultComputeHessian", flg, &flg, NULL));
509:   if (flg) {
510:     Mat H;

512:     PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &H));
513:     PetscCall(MatSetType(H, MATAIJ));
514:     PetscCall(TaoSetHessian(tao, H, H, TaoDefaultComputeHessian, NULL));
515:     PetscCall(MatDestroy(&H));
516:   }
517:   flg = PETSC_FALSE;
518:   PetscCall(PetscOptionsBool("-tao_mf_hessian", "compute matrix-free Hessian using finite differences", "TaoDefaultComputeHessianMFFD", flg, &flg, NULL));
519:   if (flg) {
520:     Mat H;

522:     PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &H));
523:     PetscCall(TaoSetHessian(tao, H, H, TaoDefaultComputeHessianMFFD, NULL));
524:     PetscCall(MatDestroy(&H));
525:   }
526:   PetscCall(PetscOptionsBool("-tao_recycle_history", "enable recycling/re-using information from the previous TaoSolve() call for some algorithms", "TaoSetRecycleHistory", flg, &flg, &found));
527:   if (found) PetscCall(TaoSetRecycleHistory(tao, flg));
528:   PetscCall(PetscOptionsEnum("-tao_subset_type", "subset type", "", TaoSubSetTypes, (PetscEnum)tao->subset_type, (PetscEnum *)&tao->subset_type, NULL));

530:   if (tao->ksp) {
531:     PetscCall(PetscOptionsBool("-tao_ksp_ew", "Use Eisentat-Walker linear system convergence test", "TaoKSPSetUseEW", tao->ksp_ewconv, &tao->ksp_ewconv, NULL));
532:     PetscCall(TaoKSPSetUseEW(tao, tao->ksp_ewconv));
533:   }

535:   PetscTryTypeMethod(tao, setfromoptions, PetscOptionsObject);

537:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
538:   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)tao, PetscOptionsObject));
539:   PetscOptionsEnd();

541:   if (tao->linesearch) PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
542:   PetscFunctionReturn(PETSC_SUCCESS);
543: }

545: /*@C
546:   TaoViewFromOptions - View a `Tao` object based on values in the options database

548:   Collective

550:   Input Parameters:
551: + A    - the  `Tao` context
552: . obj  - Optional object that provides the prefix for the options database
553: - name - command line option

555:   Level: intermediate

557: .seealso: [](ch_tao), `Tao`, `TaoView`, `PetscObjectViewFromOptions()`, `TaoCreate()`
558: @*/
559: PetscErrorCode TaoViewFromOptions(Tao A, PetscObject obj, const char name[])
560: {
561:   PetscFunctionBegin;
563:   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
564:   PetscFunctionReturn(PETSC_SUCCESS);
565: }

567: /*@C
568:   TaoView - Prints information about the `Tao` object

570:   Collective

572:   Input Parameters:
573: + tao    - the `Tao` context
574: - viewer - visualization context

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

579:   Level: beginner

581:   Notes:
582:   The available visualization contexts include
583: +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
584: -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
585:   output where only the first processor opens
586:   the file.  All other processors send their
587:   data to the first processor to print.

589: .seealso: [](ch_tao), `Tao`, `PetscViewerASCIIOpen()`
590: @*/
591: PetscErrorCode TaoView(Tao tao, PetscViewer viewer)
592: {
593:   PetscBool isascii, isstring;
594:   TaoType   type;

596:   PetscFunctionBegin;
598:   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(((PetscObject)tao)->comm, &viewer));
600:   PetscCheckSameComm(tao, 1, viewer, 2);

602:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
603:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
604:   if (isascii) {
605:     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)tao, viewer));

607:     if (tao->ops->view) {
608:       PetscCall(PetscViewerASCIIPushTab(viewer));
609:       PetscUseTypeMethod(tao, view, viewer);
610:       PetscCall(PetscViewerASCIIPopTab(viewer));
611:     }
612:     if (tao->linesearch) {
613:       PetscCall(PetscViewerASCIIPushTab(viewer));
614:       PetscCall(TaoLineSearchView(tao->linesearch, viewer));
615:       PetscCall(PetscViewerASCIIPopTab(viewer));
616:     }
617:     if (tao->ksp) {
618:       PetscCall(PetscViewerASCIIPushTab(viewer));
619:       PetscCall(KSPView(tao->ksp, viewer));
620:       PetscCall(PetscViewerASCIIPrintf(viewer, "total KSP iterations: %" PetscInt_FMT "\n", tao->ksp_tot_its));
621:       PetscCall(PetscViewerASCIIPopTab(viewer));
622:     }

624:     PetscCall(PetscViewerASCIIPushTab(viewer));

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

628:     PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: gatol=%g,", (double)tao->gatol));
629:     PetscCall(PetscViewerASCIIPrintf(viewer, " steptol=%g,", (double)tao->steptol));
630:     PetscCall(PetscViewerASCIIPrintf(viewer, " gttol=%g\n", (double)tao->gttol));
631:     PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Function/Gradient:=%g\n", (double)tao->residual));

633:     if (tao->constrained) {
634:       PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances:"));
635:       PetscCall(PetscViewerASCIIPrintf(viewer, " catol=%g,", (double)tao->catol));
636:       PetscCall(PetscViewerASCIIPrintf(viewer, " crtol=%g\n", (double)tao->crtol));
637:       PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Constraints:=%g\n", (double)tao->cnorm));
638:     }

640:     if (tao->trust < tao->steptol) {
641:       PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: steptol=%g\n", (double)tao->steptol));
642:       PetscCall(PetscViewerASCIIPrintf(viewer, "Final trust region radius:=%g\n", (double)tao->trust));
643:     }

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

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

651:     if (tao->nfuncs > 0) {
652:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function evaluations=%" PetscInt_FMT ",", tao->nfuncs));
653:       PetscCall(PetscViewerASCIIPrintf(viewer, "                max: %" PetscInt_FMT "\n", tao->max_funcs));
654:     }
655:     if (tao->ngrads > 0) {
656:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of gradient evaluations=%" PetscInt_FMT ",", tao->ngrads));
657:       PetscCall(PetscViewerASCIIPrintf(viewer, "                max: %" PetscInt_FMT "\n", tao->max_funcs));
658:     }
659:     if (tao->nfuncgrads > 0) {
660:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function/gradient evaluations=%" PetscInt_FMT ",", tao->nfuncgrads));
661:       PetscCall(PetscViewerASCIIPrintf(viewer, "    (max: %" PetscInt_FMT ")\n", tao->max_funcs));
662:     }
663:     if (tao->nhess > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Hessian evaluations=%" PetscInt_FMT "\n", tao->nhess));
664:     if (tao->nconstraints > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of constraint function evaluations=%" PetscInt_FMT "\n", tao->nconstraints));
665:     if (tao->njac > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Jacobian evaluations=%" PetscInt_FMT "\n", tao->njac));

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

726: /*@
727:   TaoSetRecycleHistory - Sets the boolean flag to enable/disable re-using
728:   iterate information from the previous `TaoSolve()`. This feature is disabled by
729:   default.

731:   Logically Collective

733:   Input Parameters:
734: + tao     - the `Tao` context
735: - recycle - boolean flag

737:   Options Database Key:
738: . -tao_recycle_history <true,false> - reuse the history

740:   Level: intermediate

742:   Notes:
743:   For conjugate gradient methods (`TAOBNCG`), this re-uses the latest search direction
744:   from the previous `TaoSolve()` call when computing the first search direction in a
745:   new solution. By default, CG methods set the first search direction to the
746:   negative gradient.

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

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

755: .seealso: [](ch_tao), `Tao`, `TaoGetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
756: @*/
757: PetscErrorCode TaoSetRecycleHistory(Tao tao, PetscBool recycle)
758: {
759:   PetscFunctionBegin;
762:   tao->recycle = recycle;
763:   PetscFunctionReturn(PETSC_SUCCESS);
764: }

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

770:   Logically Collective

772:   Input Parameter:
773: . tao - the `Tao` context

775:   Output Parameter:
776: . recycle - boolean flag

778:   Level: intermediate

780: .seealso: [](ch_tao), `Tao`, `TaoSetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
781: @*/
782: PetscErrorCode TaoGetRecycleHistory(Tao tao, PetscBool *recycle)
783: {
784:   PetscFunctionBegin;
786:   PetscAssertPointer(recycle, 2);
787:   *recycle = tao->recycle;
788:   PetscFunctionReturn(PETSC_SUCCESS);
789: }

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

794:   Logically Collective

796:   Input Parameters:
797: + tao   - the `Tao` context
798: . gatol - stop if norm of gradient is less than this
799: . grtol - stop if relative norm of gradient is less than this
800: - gttol - stop if norm of gradient is reduced by this factor

802:   Options Database Keys:
803: + -tao_gatol <gatol> - Sets gatol
804: . -tao_grtol <grtol> - Sets grtol
805: - -tao_gttol <gttol> - Sets gttol

807:   Stopping Criteria\:
808: .vb
809:   ||g(X)||                            <= gatol
810:   ||g(X)|| / |f(X)|                   <= grtol
811:   ||g(X)|| / ||g(X0)||                <= gttol
812: .ve

814:   Level: beginner

816:   Note:
817:   Use `PETSC_DEFAULT` to leave one or more tolerances unchanged.

819: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`
820: @*/
821: PetscErrorCode TaoSetTolerances(Tao tao, PetscReal gatol, PetscReal grtol, PetscReal gttol)
822: {
823:   PetscFunctionBegin;

829:   if (gatol != (PetscReal)PETSC_DEFAULT) {
830:     if (gatol < 0) {
831:       PetscCall(PetscInfo(tao, "Tried to set negative gatol -- ignored.\n"));
832:     } else {
833:       tao->gatol         = PetscMax(0, gatol);
834:       tao->gatol_changed = PETSC_TRUE;
835:     }
836:   }

838:   if (grtol != (PetscReal)PETSC_DEFAULT) {
839:     if (grtol < 0) {
840:       PetscCall(PetscInfo(tao, "Tried to set negative grtol -- ignored.\n"));
841:     } else {
842:       tao->grtol         = PetscMax(0, grtol);
843:       tao->grtol_changed = PETSC_TRUE;
844:     }
845:   }

847:   if (gttol != (PetscReal)PETSC_DEFAULT) {
848:     if (gttol < 0) {
849:       PetscCall(PetscInfo(tao, "Tried to set negative gttol -- ignored.\n"));
850:     } else {
851:       tao->gttol         = PetscMax(0, gttol);
852:       tao->gttol_changed = PETSC_TRUE;
853:     }
854:   }
855:   PetscFunctionReturn(PETSC_SUCCESS);
856: }

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

861:   Logically Collective

863:   Input Parameters:
864: + tao   - the `Tao` context
865: . catol - absolute constraint tolerance, constraint norm must be less than `catol` for used for gatol convergence criteria
866: - crtol - relative constraint tolerance, constraint norm must be less than `crtol` for used for gatol, gttol convergence criteria

868:   Options Database Keys:
869: + -tao_catol <catol> - Sets catol
870: - -tao_crtol <crtol> - Sets crtol

872:   Level: intermediate

874:   Notes:
875:   Use `PETSC_DEFAULT` to leave any tolerance unchanged.

877: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`, `TaoGetConstraintTolerances()`, `TaoSetTolerances()`
878: @*/
879: PetscErrorCode TaoSetConstraintTolerances(Tao tao, PetscReal catol, PetscReal crtol)
880: {
881:   PetscFunctionBegin;

886:   if (catol != (PetscReal)PETSC_DEFAULT) {
887:     if (catol < 0) {
888:       PetscCall(PetscInfo(tao, "Tried to set negative catol -- ignored.\n"));
889:     } else {
890:       tao->catol         = PetscMax(0, catol);
891:       tao->catol_changed = PETSC_TRUE;
892:     }
893:   }

895:   if (crtol != (PetscReal)PETSC_DEFAULT) {
896:     if (crtol < 0) {
897:       PetscCall(PetscInfo(tao, "Tried to set negative crtol -- ignored.\n"));
898:     } else {
899:       tao->crtol         = PetscMax(0, crtol);
900:       tao->crtol_changed = PETSC_TRUE;
901:     }
902:   }
903:   PetscFunctionReturn(PETSC_SUCCESS);
904: }

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

909:   Not Collective

911:   Input Parameter:
912: . tao - the `Tao` context

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

918:   Level: intermediate

920: .seealso: [](ch_tao), `Tao`, `TaoConvergedReasons`,`TaoGetTolerances()`, `TaoSetTolerances()`, `TaoSetConstraintTolerances()`
921: @*/
922: PetscErrorCode TaoGetConstraintTolerances(Tao tao, PetscReal *catol, PetscReal *crtol)
923: {
924:   PetscFunctionBegin;
926:   if (catol) *catol = tao->catol;
927:   if (crtol) *crtol = tao->crtol;
928:   PetscFunctionReturn(PETSC_SUCCESS);
929: }

931: /*@
932:   TaoSetFunctionLowerBound - Sets a bound on the solution objective value.
933:   When an approximate solution with an objective value below this number
934:   has been found, the solver will terminate.

936:   Logically Collective

938:   Input Parameters:
939: + tao  - the Tao solver context
940: - fmin - the tolerance

942:   Options Database Key:
943: . -tao_fmin <fmin> - sets the minimum function value

945:   Level: intermediate

947: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetTolerances()`
948: @*/
949: PetscErrorCode TaoSetFunctionLowerBound(Tao tao, PetscReal fmin)
950: {
951:   PetscFunctionBegin;
954:   tao->fmin         = fmin;
955:   tao->fmin_changed = PETSC_TRUE;
956:   PetscFunctionReturn(PETSC_SUCCESS);
957: }

959: /*@
960:   TaoGetFunctionLowerBound - Gets the bound on the solution objective value.
961:   When an approximate solution with an objective value below this number
962:   has been found, the solver will terminate.

964:   Not Collective

966:   Input Parameter:
967: . tao - the `Tao` solver context

969:   Output Parameter:
970: . fmin - the minimum function value

972:   Level: intermediate

974: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetFunctionLowerBound()`
975: @*/
976: PetscErrorCode TaoGetFunctionLowerBound(Tao tao, PetscReal *fmin)
977: {
978:   PetscFunctionBegin;
980:   PetscAssertPointer(fmin, 2);
981:   *fmin = tao->fmin;
982:   PetscFunctionReturn(PETSC_SUCCESS);
983: }

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

988:   Logically Collective

990:   Input Parameters:
991: + tao  - the `Tao` solver context
992: - nfcn - the maximum number of function evaluations (>=0)

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

997:   Level: intermediate

999: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoSetMaximumIterations()`
1000: @*/
1001: PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao tao, PetscInt nfcn)
1002: {
1003:   PetscFunctionBegin;
1006:   if (nfcn >= 0) {
1007:     tao->max_funcs = PetscMax(0, nfcn);
1008:   } else {
1009:     tao->max_funcs = -1;
1010:   }
1011:   tao->max_funcs_changed = PETSC_TRUE;
1012:   PetscFunctionReturn(PETSC_SUCCESS);
1013: }

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

1018:   Logically Collective

1020:   Input Parameter:
1021: . tao - the `Tao` solver context

1023:   Output Parameter:
1024: . nfcn - the maximum number of function evaluations

1026:   Level: intermediate

1028: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1029: @*/
1030: PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao tao, PetscInt *nfcn)
1031: {
1032:   PetscFunctionBegin;
1034:   PetscAssertPointer(nfcn, 2);
1035:   *nfcn = tao->max_funcs;
1036:   PetscFunctionReturn(PETSC_SUCCESS);
1037: }

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

1042:   Not Collective

1044:   Input Parameter:
1045: . tao - the `Tao` solver context

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

1050:   Level: intermediate

1052: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1053: @*/
1054: PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao tao, PetscInt *nfuncs)
1055: {
1056:   PetscFunctionBegin;
1058:   PetscAssertPointer(nfuncs, 2);
1059:   *nfuncs = PetscMax(tao->nfuncs, tao->nfuncgrads);
1060:   PetscFunctionReturn(PETSC_SUCCESS);
1061: }

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

1066:   Logically Collective

1068:   Input Parameters:
1069: + tao    - the `Tao` solver context
1070: - maxits - the maximum number of iterates (>=0)

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

1075:   Level: intermediate

1077: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoSetMaximumFunctionEvaluations()`
1078: @*/
1079: PetscErrorCode TaoSetMaximumIterations(Tao tao, PetscInt maxits)
1080: {
1081:   PetscFunctionBegin;
1084:   tao->max_it         = PetscMax(0, maxits);
1085:   tao->max_it_changed = PETSC_TRUE;
1086:   PetscFunctionReturn(PETSC_SUCCESS);
1087: }

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

1092:   Not Collective

1094:   Input Parameter:
1095: . tao - the `Tao` solver context

1097:   Output Parameter:
1098: . maxits - the maximum number of iterates

1100:   Level: intermediate

1102: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumIterations()`, `TaoGetMaximumFunctionEvaluations()`
1103: @*/
1104: PetscErrorCode TaoGetMaximumIterations(Tao tao, PetscInt *maxits)
1105: {
1106:   PetscFunctionBegin;
1108:   PetscAssertPointer(maxits, 2);
1109:   *maxits = tao->max_it;
1110:   PetscFunctionReturn(PETSC_SUCCESS);
1111: }

1113: /*@
1114:   TaoSetInitialTrustRegionRadius - Sets the initial trust region radius.

1116:   Logically Collective

1118:   Input Parameters:
1119: + tao    - a `Tao` optimization solver
1120: - radius - the trust region radius

1122:   Options Database Key:
1123: . -tao_trust0 <t0> - sets initial trust region radius

1125:   Level: intermediate

1127: .seealso: [](ch_tao), `Tao`, `TaoGetTrustRegionRadius()`, `TaoSetTrustRegionTolerance()`, `TAONTR`
1128: @*/
1129: PetscErrorCode TaoSetInitialTrustRegionRadius(Tao tao, PetscReal radius)
1130: {
1131:   PetscFunctionBegin;
1134:   tao->trust0         = PetscMax(0.0, radius);
1135:   tao->trust0_changed = PETSC_TRUE;
1136:   PetscFunctionReturn(PETSC_SUCCESS);
1137: }

1139: /*@
1140:   TaoGetInitialTrustRegionRadius - Gets the initial trust region radius.

1142:   Not Collective

1144:   Input Parameter:
1145: . tao - a `Tao` optimization solver

1147:   Output Parameter:
1148: . radius - the trust region radius

1150:   Level: intermediate

1152: .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetCurrentTrustRegionRadius()`, `TAONTR`
1153: @*/
1154: PetscErrorCode TaoGetInitialTrustRegionRadius(Tao tao, PetscReal *radius)
1155: {
1156:   PetscFunctionBegin;
1158:   PetscAssertPointer(radius, 2);
1159:   *radius = tao->trust0;
1160:   PetscFunctionReturn(PETSC_SUCCESS);
1161: }

1163: /*@
1164:   TaoGetCurrentTrustRegionRadius - Gets the current trust region radius.

1166:   Not Collective

1168:   Input Parameter:
1169: . tao - a `Tao` optimization solver

1171:   Output Parameter:
1172: . radius - the trust region radius

1174:   Level: intermediate

1176: .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetInitialTrustRegionRadius()`, `TAONTR`
1177: @*/
1178: PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao tao, PetscReal *radius)
1179: {
1180:   PetscFunctionBegin;
1182:   PetscAssertPointer(radius, 2);
1183:   *radius = tao->trust;
1184:   PetscFunctionReturn(PETSC_SUCCESS);
1185: }

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

1190:   Not Collective

1192:   Input Parameter:
1193: . tao - the `Tao` context

1195:   Output Parameters:
1196: + gatol - stop if norm of gradient is less than this
1197: . grtol - stop if relative norm of gradient is less than this
1198: - gttol - stop if norm of gradient is reduced by a this factor

1200:   Level: intermediate

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

1205: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`
1206: @*/
1207: PetscErrorCode TaoGetTolerances(Tao tao, PetscReal *gatol, PetscReal *grtol, PetscReal *gttol)
1208: {
1209:   PetscFunctionBegin;
1211:   if (gatol) *gatol = tao->gatol;
1212:   if (grtol) *grtol = tao->grtol;
1213:   if (gttol) *gttol = tao->gttol;
1214:   PetscFunctionReturn(PETSC_SUCCESS);
1215: }

1217: /*@
1218:   TaoGetKSP - Gets the linear solver used by the optimization solver.

1220:   Not Collective

1222:   Input Parameter:
1223: . tao - the `Tao` solver

1225:   Output Parameter:
1226: . ksp - the `KSP` linear solver used in the optimization solver

1228:   Level: intermediate

1230: .seealso: [](ch_tao), `Tao`, `KSP`
1231: @*/
1232: PetscErrorCode TaoGetKSP(Tao tao, KSP *ksp)
1233: {
1234:   PetscFunctionBegin;
1236:   PetscAssertPointer(ksp, 2);
1237:   *ksp = tao->ksp;
1238:   PetscFunctionReturn(PETSC_SUCCESS);
1239: }

1241: /*@
1242:   TaoGetLinearSolveIterations - Gets the total number of linear iterations
1243:   used by the `Tao` solver

1245:   Not Collective

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

1250:   Output Parameter:
1251: . lits - number of linear iterations

1253:   Level: intermediate

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

1258: .seealso: [](ch_tao), `Tao`, `TaoGetKSP()`
1259: @*/
1260: PetscErrorCode TaoGetLinearSolveIterations(Tao tao, PetscInt *lits)
1261: {
1262:   PetscFunctionBegin;
1264:   PetscAssertPointer(lits, 2);
1265:   *lits = tao->ksp_tot_its;
1266:   PetscFunctionReturn(PETSC_SUCCESS);
1267: }

1269: /*@
1270:   TaoGetLineSearch - Gets the line search used by the optimization solver.

1272:   Not Collective

1274:   Input Parameter:
1275: . tao - the `Tao` solver

1277:   Output Parameter:
1278: . ls - the line search used in the optimization solver

1280:   Level: intermediate

1282: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`
1283: @*/
1284: PetscErrorCode TaoGetLineSearch(Tao tao, TaoLineSearch *ls)
1285: {
1286:   PetscFunctionBegin;
1288:   PetscAssertPointer(ls, 2);
1289:   *ls = tao->linesearch;
1290:   PetscFunctionReturn(PETSC_SUCCESS);
1291: }

1293: /*@
1294:   TaoAddLineSearchCounts - Adds the number of function evaluations spent
1295:   in the line search to the running total.

1297:   Input Parameters:
1298: . tao - the `Tao` solver

1300:   Level: developer

1302: .seealso: [](ch_tao), `Tao`, `TaoGetLineSearch()`, `TaoLineSearchApply()`
1303: @*/
1304: PetscErrorCode TaoAddLineSearchCounts(Tao tao)
1305: {
1306:   PetscBool flg;
1307:   PetscInt  nfeval, ngeval, nfgeval;

1309:   PetscFunctionBegin;
1311:   if (tao->linesearch) {
1312:     PetscCall(TaoLineSearchIsUsingTaoRoutines(tao->linesearch, &flg));
1313:     if (!flg) {
1314:       PetscCall(TaoLineSearchGetNumberFunctionEvaluations(tao->linesearch, &nfeval, &ngeval, &nfgeval));
1315:       tao->nfuncs += nfeval;
1316:       tao->ngrads += ngeval;
1317:       tao->nfuncgrads += nfgeval;
1318:     }
1319:   }
1320:   PetscFunctionReturn(PETSC_SUCCESS);
1321: }

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

1326:   Not Collective

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

1331:   Output Parameter:
1332: . X - the current solution

1334:   Level: intermediate

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

1339: .seealso: [](ch_tao), `Tao`, `TaoSetSolution()`, `TaoSolve()`
1340: @*/
1341: PetscErrorCode TaoGetSolution(Tao tao, Vec *X)
1342: {
1343:   PetscFunctionBegin;
1345:   PetscAssertPointer(X, 2);
1346:   *X = tao->solution;
1347:   PetscFunctionReturn(PETSC_SUCCESS);
1348: }

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

1355:   Collective

1357:   Input Parameter:
1358: . tao - the `Tao` context

1360:   Level: developer

1362: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
1363: @*/
1364: PetscErrorCode TaoResetStatistics(Tao tao)
1365: {
1366:   PetscFunctionBegin;
1368:   tao->niter        = 0;
1369:   tao->nfuncs       = 0;
1370:   tao->nfuncgrads   = 0;
1371:   tao->ngrads       = 0;
1372:   tao->nhess        = 0;
1373:   tao->njac         = 0;
1374:   tao->nconstraints = 0;
1375:   tao->ksp_its      = 0;
1376:   tao->ksp_tot_its  = 0;
1377:   tao->reason       = TAO_CONTINUE_ITERATING;
1378:   tao->residual     = 0.0;
1379:   tao->cnorm        = 0.0;
1380:   tao->step         = 0.0;
1381:   tao->lsflag       = PETSC_FALSE;
1382:   if (tao->hist_reset) tao->hist_len = 0;
1383:   PetscFunctionReturn(PETSC_SUCCESS);
1384: }

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

1391:   Logically Collective

1393:   Input Parameters:
1394: + tao  - The `Tao` solver context
1395: - func - The function

1397:   Calling sequence of `func`:
1398: + tao - the optimizer context
1399: - ctx - The current step of the iteration

1401:   Level: advanced

1403: .seealso: [](ch_tao), `Tao`, `TaoSolve()`
1404: @*/
1405: PetscErrorCode TaoSetUpdate(Tao tao, PetscErrorCode (*func)(Tao, PetscInt, void *), void *ctx)
1406: {
1407:   PetscFunctionBegin;
1409:   tao->ops->update = func;
1410:   tao->user_update = ctx;
1411:   PetscFunctionReturn(PETSC_SUCCESS);
1412: }

1414: /*@C
1415:   TaoSetConvergenceTest - Sets the function that is to be used to test
1416:   for convergence o fthe iterative minimization solution.  The new convergence
1417:   testing routine will replace Tao's default convergence test.

1419:   Logically Collective

1421:   Input Parameters:
1422: + tao  - the `Tao` object
1423: . conv - the routine to test for convergence
1424: - ctx  - [optional] context for private data for the convergence routine
1425:         (may be `NULL`)

1427:   Calling sequence of `conv`:
1428: + tao - the `Tao` object
1429: - ctx - [optional] convergence context

1431:   Level: advanced

1433:   Note:
1434:   The new convergence testing routine should call `TaoSetConvergedReason()`.

1436: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetConvergedReason()`, `TaoGetSolutionStatus()`, `TaoGetTolerances()`, `TaoMonitorSet()`
1437: @*/
1438: PetscErrorCode TaoSetConvergenceTest(Tao tao, PetscErrorCode (*conv)(Tao, void *), void *ctx)
1439: {
1440:   PetscFunctionBegin;
1442:   tao->ops->convergencetest = conv;
1443:   tao->cnvP                 = ctx;
1444:   PetscFunctionReturn(PETSC_SUCCESS);
1445: }

1447: /*@C
1448:   TaoMonitorSet - Sets an additional function that is to be used at every
1449:   iteration of the solver to display the iteration's
1450:   progress.

1452:   Logically Collective

1454:   Input Parameters:
1455: + tao  - the `Tao` solver context
1456: . func - monitoring routine
1457: . ctx  - [optional] user-defined context for private data for the monitor routine (may be `NULL`)
1458: - dest - [optional] function to destroy the context when the `Tao` is destroyed

1460:   Calling sequence of `func`:
1461: + tao - the `Tao` solver context
1462: - ctx - [optional] monitoring context

1464:   Calling sequence of `dest`:
1465: . ctx - monitoring context

1467:   Level: intermediate

1469:   Notes:
1470:   See `TaoSetFromOptions()` for a monitoring options.

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

1476:   Fortran Notes:
1477:   Only one monitor function may be set

1479: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoMonitorDefault()`, `TaoMonitorCancel()`, `TaoSetDestroyRoutine()`, `TaoView()`
1480: @*/
1481: PetscErrorCode TaoMonitorSet(Tao tao, PetscErrorCode (*func)(Tao, void *), void *ctx, PetscErrorCode (*dest)(void **))
1482: {
1483:   PetscInt  i;
1484:   PetscBool identical;

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

1490:   for (i = 0; i < tao->numbermonitors; i++) {
1491:     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))func, ctx, dest, (PetscErrorCode(*)(void))tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i], &identical));
1492:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
1493:   }
1494:   tao->monitor[tao->numbermonitors]        = func;
1495:   tao->monitorcontext[tao->numbermonitors] = (void *)ctx;
1496:   tao->monitordestroy[tao->numbermonitors] = dest;
1497:   ++tao->numbermonitors;
1498:   PetscFunctionReturn(PETSC_SUCCESS);
1499: }

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

1504:   Logically Collective

1506:   Input Parameter:
1507: . tao - the `Tao` solver context

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

1514:   Level: advanced

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

1519: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1520: @*/
1521: PetscErrorCode TaoMonitorCancel(Tao tao)
1522: {
1523:   PetscInt i;

1525:   PetscFunctionBegin;
1527:   for (i = 0; i < tao->numbermonitors; i++) {
1528:     if (tao->monitordestroy[i]) PetscCall((*tao->monitordestroy[i])(&tao->monitorcontext[i]));
1529:   }
1530:   tao->numbermonitors = 0;
1531:   PetscFunctionReturn(PETSC_SUCCESS);
1532: }

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

1537:   Collective

1539:   Input Parameters:
1540: + tao - the `Tao` context
1541: - ctx - `PetscViewer` context or `NULL`

1543:   Options Database Key:
1544: . -tao_monitor - turn on default monitoring

1546:   Level: advanced

1548:   Note:
1549:   This monitor prints the function value and gradient
1550:   norm at each iteration.

1552: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1553: @*/
1554: PetscErrorCode TaoMonitorDefault(Tao tao, void *ctx)
1555: {
1556:   PetscInt    its, tabs;
1557:   PetscReal   fct, gnorm;
1558:   PetscViewer viewer = (PetscViewer)ctx;

1560:   PetscFunctionBegin;
1563:   its   = tao->niter;
1564:   fct   = tao->fc;
1565:   gnorm = tao->residual;
1566:   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1567:   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1568:   if (its == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
1569:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Iteration information for %s solve.\n", ((PetscObject)tao)->prefix));
1570:     tao->header_printed = PETSC_TRUE;
1571:   }
1572:   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " TAO,", its));
1573:   PetscCall(PetscViewerASCIIPrintf(viewer, "  Function value: %g,", (double)fct));
1574:   if (gnorm >= PETSC_INFINITY) {
1575:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: Inf \n"));
1576:   } else {
1577:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g \n", (double)gnorm));
1578:   }
1579:   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1580:   PetscFunctionReturn(PETSC_SUCCESS);
1581: }

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

1586:   Collective

1588:   Input Parameters:
1589: + tao - the `Tao` context
1590: - ctx - `PetscViewer` context or `NULL`

1592:   Options Database Key:
1593: . -tao_monitor_globalization - turn on monitoring with globalization information

1595:   Level: advanced

1597:   Note:
1598:   This monitor prints the function value and gradient norm at each
1599:   iteration, as well as the step size and trust radius. Note that the
1600:   step size and trust radius may be the same for some algorithms.

1602: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1603: @*/
1604: PetscErrorCode TaoMonitorGlobalization(Tao tao, void *ctx)
1605: {
1606:   PetscInt    its, tabs;
1607:   PetscReal   fct, gnorm, stp, tr;
1608:   PetscViewer viewer = (PetscViewer)ctx;

1610:   PetscFunctionBegin;
1613:   its   = tao->niter;
1614:   fct   = tao->fc;
1615:   gnorm = tao->residual;
1616:   stp   = tao->step;
1617:   tr    = tao->trust;
1618:   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1619:   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1620:   if (its == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
1621:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Iteration information for %s solve.\n", ((PetscObject)tao)->prefix));
1622:     tao->header_printed = PETSC_TRUE;
1623:   }
1624:   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " TAO,", its));
1625:   PetscCall(PetscViewerASCIIPrintf(viewer, "  Function value: %g,", (double)fct));
1626:   if (gnorm >= PETSC_INFINITY) {
1627:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: Inf,"));
1628:   } else {
1629:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g,", (double)gnorm));
1630:   }
1631:   PetscCall(PetscViewerASCIIPrintf(viewer, "  Step: %g,  Trust: %g\n", (double)stp, (double)tr));
1632:   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1633:   PetscFunctionReturn(PETSC_SUCCESS);
1634: }

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

1639:   Collective

1641:   Input Parameters:
1642: + tao - the `Tao` context
1643: - ctx - `PetscViewer` context of type `PETSCVIEWERASCII`

1645:   Options Database Key:
1646: . -tao_monitor_short - turn on default short monitoring

1648:   Level: advanced

1650:   Note:
1651:   Same as `TaoMonitorDefault()` except
1652:   it prints fewer digits of the residual as the residual gets smaller.
1653:   This is because the later digits are meaningless and are often
1654:   different on different machines; by using this routine different
1655:   machines will usually generate the same output.

1657: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1658: @*/
1659: PetscErrorCode TaoMonitorDefaultShort(Tao tao, void *ctx)
1660: {
1661:   PetscInt    its, tabs;
1662:   PetscReal   fct, gnorm;
1663:   PetscViewer viewer = (PetscViewer)ctx;

1665:   PetscFunctionBegin;
1668:   its   = tao->niter;
1669:   fct   = tao->fc;
1670:   gnorm = tao->residual;
1671:   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1672:   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1673:   PetscCall(PetscViewerASCIIPrintf(viewer, "iter = %3" PetscInt_FMT ",", its));
1674:   PetscCall(PetscViewerASCIIPrintf(viewer, " Function value %g,", (double)fct));
1675:   if (gnorm >= PETSC_INFINITY) {
1676:     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: Inf \n"));
1677:   } else if (gnorm > 1.e-6) {
1678:     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: %g \n", (double)gnorm));
1679:   } else if (gnorm > 1.e-11) {
1680:     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: < 1.0e-6 \n"));
1681:   } else {
1682:     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: < 1.0e-11 \n"));
1683:   }
1684:   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1685:   PetscFunctionReturn(PETSC_SUCCESS);
1686: }

1688: /*@
1689:   TaoMonitorConstraintNorm - same as `TaoMonitorDefault()` except
1690:   it prints the norm of the constraint function.

1692:   Collective

1694:   Input Parameters:
1695: + tao - the `Tao` context
1696: - ctx - `PetscViewer` context or `NULL`

1698:   Options Database Key:
1699: . -tao_monitor_constraint_norm - monitor the constraints

1701:   Level: advanced

1703: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1704: @*/
1705: PetscErrorCode TaoMonitorConstraintNorm(Tao tao, void *ctx)
1706: {
1707:   PetscInt    its, tabs;
1708:   PetscReal   fct, gnorm;
1709:   PetscViewer viewer = (PetscViewer)ctx;

1711:   PetscFunctionBegin;
1714:   its   = tao->niter;
1715:   fct   = tao->fc;
1716:   gnorm = tao->residual;
1717:   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1718:   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1719:   PetscCall(PetscViewerASCIIPrintf(viewer, "iter = %" PetscInt_FMT ",", its));
1720:   PetscCall(PetscViewerASCIIPrintf(viewer, " Function value: %g,", (double)fct));
1721:   PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g ", (double)gnorm));
1722:   PetscCall(PetscViewerASCIIPrintf(viewer, "  Constraint: %g \n", (double)tao->cnorm));
1723:   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1724:   PetscFunctionReturn(PETSC_SUCCESS);
1725: }

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

1730:   Collective

1732:   Input Parameters:
1733: + tao - the `Tao` context
1734: - ctx - `PetscViewer` context or `NULL`

1736:   Options Database Key:
1737: . -tao_monitor_solution - view the solution

1739:   Level: advanced

1741: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1742: @*/
1743: PetscErrorCode TaoMonitorSolution(Tao tao, void *ctx)
1744: {
1745:   PetscViewer viewer = (PetscViewer)ctx;

1747:   PetscFunctionBegin;
1750:   PetscCall(VecView(tao->solution, viewer));
1751:   PetscFunctionReturn(PETSC_SUCCESS);
1752: }

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

1757:   Collective

1759:   Input Parameters:
1760: + tao - the `Tao` context
1761: - ctx - `PetscViewer` context or `NULL`

1763:   Options Database Key:
1764: . -tao_monitor_gradient - view the gradient at each iteration

1766:   Level: advanced

1768: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1769: @*/
1770: PetscErrorCode TaoMonitorGradient(Tao tao, void *ctx)
1771: {
1772:   PetscViewer viewer = (PetscViewer)ctx;

1774:   PetscFunctionBegin;
1777:   PetscCall(VecView(tao->gradient, viewer));
1778:   PetscFunctionReturn(PETSC_SUCCESS);
1779: }

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

1784:   Collective

1786:   Input Parameters:
1787: + tao - the `Tao` context
1788: - ctx - `PetscViewer` context or `NULL`

1790:   Options Database Key:
1791: . -tao_monitor_step - view the step vector at each iteration

1793:   Level: advanced

1795: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1796: @*/
1797: PetscErrorCode TaoMonitorStep(Tao tao, void *ctx)
1798: {
1799:   PetscViewer viewer = (PetscViewer)ctx;

1801:   PetscFunctionBegin;
1804:   PetscCall(VecView(tao->stepdirection, viewer));
1805:   PetscFunctionReturn(PETSC_SUCCESS);
1806: }

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

1811:   Collective

1813:   Input Parameters:
1814: + tao - the `Tao` context
1815: - ctx - `TaoMonitorDraw` context

1817:   Options Database Key:
1818: . -tao_monitor_solution_draw - draw the solution at each iteration

1820:   Level: advanced

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

1826: .seealso: [](ch_tao), `Tao`, `TaoMonitorSolution()`, `TaoMonitorSet()`, `TaoMonitorGradientDraw()`, `TaoMonitorDrawCtxCreate()`,
1827:           `TaoMonitorDrawCtxDestroy()`
1828: @*/
1829: PetscErrorCode TaoMonitorSolutionDraw(Tao tao, void *ctx)
1830: {
1831:   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;

1833:   PetscFunctionBegin;
1835:   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(PETSC_SUCCESS);
1836:   PetscCall(VecView(tao->solution, ictx->viewer));
1837:   PetscFunctionReturn(PETSC_SUCCESS);
1838: }

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

1843:   Collective

1845:   Input Parameters:
1846: + tao - the `Tao` context
1847: - ctx - `PetscViewer` context

1849:   Options Database Key:
1850: . -tao_monitor_gradient_draw - draw the gradient at each iteration

1852:   Level: advanced

1854: .seealso: [](ch_tao), `Tao`, `TaoMonitorGradient()`, `TaoMonitorSet()`, `TaoMonitorSolutionDraw()`
1855: @*/
1856: PetscErrorCode TaoMonitorGradientDraw(Tao tao, void *ctx)
1857: {
1858:   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;

1860:   PetscFunctionBegin;
1862:   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(PETSC_SUCCESS);
1863:   PetscCall(VecView(tao->gradient, ictx->viewer));
1864:   PetscFunctionReturn(PETSC_SUCCESS);
1865: }

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

1870:   Collective

1872:   Input Parameters:
1873: + tao - the `Tao` context
1874: - ctx - the `PetscViewer` context

1876:   Options Database Key:
1877: . -tao_monitor_step_draw - draw the step direction at each iteration

1879:   Level: advanced

1881: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorSolutionDraw`
1882: @*/
1883: PetscErrorCode TaoMonitorStepDraw(Tao tao, void *ctx)
1884: {
1885:   PetscViewer viewer = (PetscViewer)ctx;

1887:   PetscFunctionBegin;
1890:   PetscCall(VecView(tao->stepdirection, viewer));
1891:   PetscFunctionReturn(PETSC_SUCCESS);
1892: }

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

1897:   Collective

1899:   Input Parameters:
1900: + tao - the `Tao` context
1901: - ctx - the `PetscViewer` context or `NULL`

1903:   Options Database Key:
1904: . -tao_monitor_ls_residual - view the residual at each iteration

1906:   Level: advanced

1908: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1909: @*/
1910: PetscErrorCode TaoMonitorResidual(Tao tao, void *ctx)
1911: {
1912:   PetscViewer viewer = (PetscViewer)ctx;

1914:   PetscFunctionBegin;
1917:   PetscCall(VecView(tao->ls_res, viewer));
1918:   PetscFunctionReturn(PETSC_SUCCESS);
1919: }

1921: /*@
1922:   TaoDefaultConvergenceTest - Determines whether the solver should continue iterating
1923:   or terminate.

1925:   Collective

1927:   Input Parameters:
1928: + tao   - the `Tao` context
1929: - dummy - unused dummy context

1931:   Level: developer

1933:   Notes:
1934:   This routine checks the residual in the optimality conditions, the
1935:   relative residual in the optimity conditions, the number of function
1936:   evaluations, and the function value to test convergence.  Some
1937:   solvers may use different convergence routines.

1939: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoGetConvergedReason()`, `TaoSetConvergedReason()`
1940: @*/
1941: PetscErrorCode TaoDefaultConvergenceTest(Tao tao, void *dummy)
1942: {
1943:   PetscInt           niter = tao->niter, nfuncs = PetscMax(tao->nfuncs, tao->nfuncgrads);
1944:   PetscInt           max_funcs = tao->max_funcs;
1945:   PetscReal          gnorm = tao->residual, gnorm0 = tao->gnorm0;
1946:   PetscReal          f = tao->fc, steptol = tao->steptol, trradius = tao->step;
1947:   PetscReal          gatol = tao->gatol, grtol = tao->grtol, gttol = tao->gttol;
1948:   PetscReal          catol = tao->catol, crtol = tao->crtol;
1949:   PetscReal          fmin = tao->fmin, cnorm = tao->cnorm;
1950:   TaoConvergedReason reason = tao->reason;

1952:   PetscFunctionBegin;
1954:   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);

1956:   if (PetscIsInfOrNanReal(f)) {
1957:     PetscCall(PetscInfo(tao, "Failed to converged, function value is Inf or NaN\n"));
1958:     reason = TAO_DIVERGED_NAN;
1959:   } else if (f <= fmin && cnorm <= catol) {
1960:     PetscCall(PetscInfo(tao, "Converged due to function value %g < minimum function value %g\n", (double)f, (double)fmin));
1961:     reason = TAO_CONVERGED_MINF;
1962:   } else if (gnorm <= gatol && cnorm <= catol) {
1963:     PetscCall(PetscInfo(tao, "Converged due to residual norm ||g(X)||=%g < %g\n", (double)gnorm, (double)gatol));
1964:     reason = TAO_CONVERGED_GATOL;
1965:   } else if (f != 0 && PetscAbsReal(gnorm / f) <= grtol && cnorm <= crtol) {
1966:     PetscCall(PetscInfo(tao, "Converged due to residual ||g(X)||/|f(X)| =%g < %g\n", (double)(gnorm / f), (double)grtol));
1967:     reason = TAO_CONVERGED_GRTOL;
1968:   } else if (gnorm0 != 0 && ((gttol == 0 && gnorm == 0) || gnorm / gnorm0 < gttol) && cnorm <= crtol) {
1969:     PetscCall(PetscInfo(tao, "Converged due to relative residual norm ||g(X)||/||g(X0)|| = %g < %g\n", (double)(gnorm / gnorm0), (double)gttol));
1970:     reason = TAO_CONVERGED_GTTOL;
1971:   } else if (max_funcs >= 0 && nfuncs > max_funcs) {
1972:     PetscCall(PetscInfo(tao, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", nfuncs, max_funcs));
1973:     reason = TAO_DIVERGED_MAXFCN;
1974:   } else if (tao->lsflag != 0) {
1975:     PetscCall(PetscInfo(tao, "Tao Line Search failure.\n"));
1976:     reason = TAO_DIVERGED_LS_FAILURE;
1977:   } else if (trradius < steptol && niter > 0) {
1978:     PetscCall(PetscInfo(tao, "Trust region/step size too small: %g < %g\n", (double)trradius, (double)steptol));
1979:     reason = TAO_CONVERGED_STEPTOL;
1980:   } else if (niter >= tao->max_it) {
1981:     PetscCall(PetscInfo(tao, "Exceeded maximum number of iterations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", niter, tao->max_it));
1982:     reason = TAO_DIVERGED_MAXITS;
1983:   } else {
1984:     reason = TAO_CONTINUE_ITERATING;
1985:   }
1986:   tao->reason = reason;
1987:   PetscFunctionReturn(PETSC_SUCCESS);
1988: }

1990: /*@C
1991:   TaoSetOptionsPrefix - Sets the prefix used for searching for all
1992:   Tao options in the database.

1994:   Logically Collective

1996:   Input Parameters:
1997: + tao - the `Tao` context
1998: - p   - the prefix string to prepend to all Tao option requests

2000:   Level: advanced

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

2006:   For example, to distinguish between the runtime options for two
2007:   different Tao solvers, one could call
2008: .vb
2009:       TaoSetOptionsPrefix(tao1,"sys1_")
2010:       TaoSetOptionsPrefix(tao2,"sys2_")
2011: .ve

2013:   This would enable use of different options for each system, such as
2014: .vb
2015:       -sys1_tao_method blmvm -sys1_tao_grtol 1.e-3
2016:       -sys2_tao_method lmvm  -sys2_tao_grtol 1.e-4
2017: .ve

2019: .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoAppendOptionsPrefix()`, `TaoGetOptionsPrefix()`
2020: @*/
2021: PetscErrorCode TaoSetOptionsPrefix(Tao tao, const char p[])
2022: {
2023:   PetscFunctionBegin;
2025:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tao, p));
2026:   if (tao->linesearch) PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, p));
2027:   if (tao->ksp) PetscCall(KSPSetOptionsPrefix(tao->ksp, p));
2028:   PetscFunctionReturn(PETSC_SUCCESS);
2029: }

2031: /*@C
2032:   TaoAppendOptionsPrefix - Appends to the prefix used for searching for all Tao options in the database.

2034:   Logically Collective

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

2040:   Level: advanced

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

2046: .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoSetOptionsPrefix()`, `TaoGetOptionsPrefix()`
2047: @*/
2048: PetscErrorCode TaoAppendOptionsPrefix(Tao tao, const char p[])
2049: {
2050:   PetscFunctionBegin;
2052:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao, p));
2053:   if (tao->linesearch) PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao->linesearch, p));
2054:   if (tao->ksp) PetscCall(KSPAppendOptionsPrefix(tao->ksp, p));
2055:   PetscFunctionReturn(PETSC_SUCCESS);
2056: }

2058: /*@C
2059:   TaoGetOptionsPrefix - Gets the prefix used for searching for all
2060:   Tao options in the database

2062:   Not Collective

2064:   Input Parameter:
2065: . tao - the `Tao` context

2067:   Output Parameter:
2068: . p - pointer to the prefix string used is returned

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

2073:   Level: advanced

2075: .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoSetOptionsPrefix()`, `TaoAppendOptionsPrefix()`
2076: @*/
2077: PetscErrorCode TaoGetOptionsPrefix(Tao tao, const char *p[])
2078: {
2079:   PetscFunctionBegin;
2081:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)tao, p));
2082:   PetscFunctionReturn(PETSC_SUCCESS);
2083: }

2085: /*@C
2086:   TaoSetType - Sets the `TaoType` for the minimization solver.

2088:   Collective

2090:   Input Parameters:
2091: + tao  - the `Tao` solver context
2092: - type - a known method

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

2098:   Level: intermediate

2100: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoGetType()`, `TaoType`
2101: @*/
2102: PetscErrorCode TaoSetType(Tao tao, TaoType type)
2103: {
2104:   PetscErrorCode (*create_xxx)(Tao);
2105:   PetscBool issame;

2107:   PetscFunctionBegin;

2110:   PetscCall(PetscObjectTypeCompare((PetscObject)tao, type, &issame));
2111:   if (issame) PetscFunctionReturn(PETSC_SUCCESS);

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

2116:   /* Destroy the existing solver information */
2117:   PetscTryTypeMethod(tao, destroy);
2118:   PetscCall(KSPDestroy(&tao->ksp));
2119:   PetscCall(TaoLineSearchDestroy(&tao->linesearch));
2120:   tao->ops->setup          = NULL;
2121:   tao->ops->solve          = NULL;
2122:   tao->ops->view           = NULL;
2123:   tao->ops->setfromoptions = NULL;
2124:   tao->ops->destroy        = NULL;

2126:   tao->setupcalled = PETSC_FALSE;

2128:   PetscCall((*create_xxx)(tao));
2129:   PetscCall(PetscObjectChangeTypeName((PetscObject)tao, type));
2130:   PetscFunctionReturn(PETSC_SUCCESS);
2131: }

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

2136:   Not Collective

2138:   Input Parameters:
2139: + sname - name of a new user-defined solver
2140: - func  - routine to Create method context

2142:   Example Usage:
2143: .vb
2144:    TaoRegister("my_solver", MySolverCreate);
2145: .ve

2147:   Then, your solver can be chosen with the procedural interface via
2148: $     TaoSetType(tao, "my_solver")
2149:   or at runtime via the option
2150: $     -tao_type my_solver

2152:   Level: advanced

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

2157: .seealso: [](ch_tao), `Tao`, `TaoSetType()`, `TaoRegisterAll()`, `TaoRegisterDestroy()`
2158: @*/
2159: PetscErrorCode TaoRegister(const char sname[], PetscErrorCode (*func)(Tao))
2160: {
2161:   PetscFunctionBegin;
2162:   PetscCall(TaoInitializePackage());
2163:   PetscCall(PetscFunctionListAdd(&TaoList, sname, (void (*)(void))func));
2164:   PetscFunctionReturn(PETSC_SUCCESS);
2165: }

2167: /*@C
2168:   TaoRegisterDestroy - Frees the list of minimization solvers that were
2169:   registered by `TaoRegister()`.

2171:   Not Collective

2173:   Level: advanced

2175: .seealso: [](ch_tao), `Tao`, `TaoRegisterAll()`, `TaoRegister()`
2176: @*/
2177: PetscErrorCode TaoRegisterDestroy(void)
2178: {
2179:   PetscFunctionBegin;
2180:   PetscCall(PetscFunctionListDestroy(&TaoList));
2181:   TaoRegisterAllCalled = PETSC_FALSE;
2182:   PetscFunctionReturn(PETSC_SUCCESS);
2183: }

2185: /*@
2186:   TaoGetIterationNumber - Gets the number of `TaoSolve()` iterations completed
2187:   at this time.

2189:   Not Collective

2191:   Input Parameter:
2192: . tao - the `Tao` context

2194:   Output Parameter:
2195: . iter - iteration number

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

2200:   Level: intermediate

2202: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`, `TaoGetResidualNorm()`, `TaoGetObjective()`
2203: @*/
2204: PetscErrorCode TaoGetIterationNumber(Tao tao, PetscInt *iter)
2205: {
2206:   PetscFunctionBegin;
2208:   PetscAssertPointer(iter, 2);
2209:   *iter = tao->niter;
2210:   PetscFunctionReturn(PETSC_SUCCESS);
2211: }

2213: /*@
2214:   TaoGetResidualNorm - Gets the current value of the norm of the residual (gradient)
2215:   at this time.

2217:   Not Collective

2219:   Input Parameter:
2220: . tao - the `Tao` context

2222:   Output Parameter:
2223: . value - the current value

2225:   Level: intermediate

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

2231: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`, `TaoGetIterationNumber()`, `TaoGetObjective()`
2232: @*/
2233: PetscErrorCode TaoGetResidualNorm(Tao tao, PetscReal *value)
2234: {
2235:   PetscFunctionBegin;
2237:   PetscAssertPointer(value, 2);
2238:   *value = tao->residual;
2239:   PetscFunctionReturn(PETSC_SUCCESS);
2240: }

2242: /*@
2243:   TaoSetIterationNumber - Sets the current iteration number.

2245:   Logically Collective

2247:   Input Parameters:
2248: + tao  - the `Tao` context
2249: - iter - iteration number

2251:   Level: developer

2253: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
2254: @*/
2255: PetscErrorCode TaoSetIterationNumber(Tao tao, PetscInt iter)
2256: {
2257:   PetscFunctionBegin;
2260:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2261:   tao->niter = iter;
2262:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2263:   PetscFunctionReturn(PETSC_SUCCESS);
2264: }

2266: /*@
2267:   TaoGetTotalIterationNumber - Gets the total number of `TaoSolve()` iterations
2268:   completed. This number keeps accumulating if multiple solves
2269:   are called with the `Tao` object.

2271:   Not Collective

2273:   Input Parameter:
2274: . tao - the `Tao` context

2276:   Output Parameter:
2277: . iter - number of iterations

2279:   Level: intermediate

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

2285: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
2286: @*/
2287: PetscErrorCode TaoGetTotalIterationNumber(Tao tao, PetscInt *iter)
2288: {
2289:   PetscFunctionBegin;
2291:   PetscAssertPointer(iter, 2);
2292:   *iter = tao->ntotalits;
2293:   PetscFunctionReturn(PETSC_SUCCESS);
2294: }

2296: /*@
2297:   TaoSetTotalIterationNumber - Sets the current total iteration number.

2299:   Logically Collective

2301:   Input Parameters:
2302: + tao  - the `Tao` context
2303: - iter - the iteration number

2305:   Level: developer

2307: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
2308: @*/
2309: PetscErrorCode TaoSetTotalIterationNumber(Tao tao, PetscInt iter)
2310: {
2311:   PetscFunctionBegin;
2314:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2315:   tao->ntotalits = iter;
2316:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2317:   PetscFunctionReturn(PETSC_SUCCESS);
2318: }

2320: /*@
2321:   TaoSetConvergedReason - Sets the termination flag on a `Tao` object

2323:   Logically Collective

2325:   Input Parameters:
2326: + tao    - the `Tao` context
2327: - reason - the `TaoConvergedReason`

2329:   Level: intermediate

2331: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`
2332: @*/
2333: PetscErrorCode TaoSetConvergedReason(Tao tao, TaoConvergedReason reason)
2334: {
2335:   PetscFunctionBegin;
2338:   tao->reason = reason;
2339:   PetscFunctionReturn(PETSC_SUCCESS);
2340: }

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

2345:   Not Collective

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

2350:   Output Parameter:
2351: . reason - value of `TaoConvergedReason`

2353:   Level: intermediate

2355: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetConvergenceTest()`, `TaoSetTolerances()`
2356: @*/
2357: PetscErrorCode TaoGetConvergedReason(Tao tao, TaoConvergedReason *reason)
2358: {
2359:   PetscFunctionBegin;
2361:   PetscAssertPointer(reason, 2);
2362:   *reason = tao->reason;
2363:   PetscFunctionReturn(PETSC_SUCCESS);
2364: }

2366: /*@
2367:   TaoGetSolutionStatus - Get the current iterate, objective value,
2368:   residual, infeasibility, and termination from a `Tao` object

2370:   Not Collective

2372:   Input Parameter:
2373: . tao - the `Tao` context

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

2383:   Level: intermediate

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

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

2390: .seealso: [](ch_tao), `TaoMonitor()`, `TaoGetConvergedReason()`
2391: @*/
2392: PetscErrorCode TaoGetSolutionStatus(Tao tao, PetscInt *its, PetscReal *f, PetscReal *gnorm, PetscReal *cnorm, PetscReal *xdiff, TaoConvergedReason *reason)
2393: {
2394:   PetscFunctionBegin;
2396:   if (its) *its = tao->niter;
2397:   if (f) *f = tao->fc;
2398:   if (gnorm) *gnorm = tao->residual;
2399:   if (cnorm) *cnorm = tao->cnorm;
2400:   if (reason) *reason = tao->reason;
2401:   if (xdiff) *xdiff = tao->step;
2402:   PetscFunctionReturn(PETSC_SUCCESS);
2403: }

2405: /*@C
2406:   TaoGetType - Gets the current `TaoType` being used in the `Tao` object

2408:   Not Collective

2410:   Input Parameter:
2411: . tao - the `Tao` solver context

2413:   Output Parameter:
2414: . type - the `TaoType`

2416:   Level: intermediate

2418: .seealso: [](ch_tao), `Tao`, `TaoType`, `TaoSetType()`
2419: @*/
2420: PetscErrorCode TaoGetType(Tao tao, TaoType *type)
2421: {
2422:   PetscFunctionBegin;
2424:   PetscAssertPointer(type, 2);
2425:   *type = ((PetscObject)tao)->type_name;
2426:   PetscFunctionReturn(PETSC_SUCCESS);
2427: }

2429: /*@C
2430:   TaoMonitor - Monitor the solver and the current solution.  This
2431:   routine will record the iteration number and residual statistics,
2432:   and call any monitors specified by the user.

2434:   Input Parameters:
2435: + tao        - the `Tao` context
2436: . its        - the current iterate number (>=0)
2437: . f          - the current objective function value
2438: . res        - the gradient norm, square root of the duality gap, or other measure indicating distance from optimality.  This measure will be recorded and
2439:           used for some termination tests.
2440: . cnorm      - the infeasibility of the current solution with regard to the constraints.
2441: - steplength - multiple of the step direction added to the previous iterate.

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

2446:   Level: developer

2448: .seealso: [](ch_tao), `Tao`, `TaoGetConvergedReason()`, `TaoMonitorDefault()`, `TaoMonitorSet()`
2449: @*/
2450: PetscErrorCode TaoMonitor(Tao tao, PetscInt its, PetscReal f, PetscReal res, PetscReal cnorm, PetscReal steplength)
2451: {
2452:   PetscInt i;

2454:   PetscFunctionBegin;
2456:   tao->fc       = f;
2457:   tao->residual = res;
2458:   tao->cnorm    = cnorm;
2459:   tao->step     = steplength;
2460:   if (!its) {
2461:     tao->cnorm0 = cnorm;
2462:     tao->gnorm0 = res;
2463:   }
2464:   PetscCall(VecLockReadPush(tao->solution));
2465:   for (i = 0; i < tao->numbermonitors; i++) PetscCall((*tao->monitor[i])(tao, tao->monitorcontext[i]));
2466:   PetscCall(VecLockReadPop(tao->solution));
2467:   PetscFunctionReturn(PETSC_SUCCESS);
2468: }

2470: /*@
2471:   TaoSetConvergenceHistory - Sets the array used to hold the convergence history.

2473:   Logically Collective

2475:   Input Parameters:
2476: + tao   - the `Tao` solver context
2477: . obj   - array to hold objective value history
2478: . resid - array to hold residual history
2479: . cnorm - array to hold constraint violation history
2480: . lits  - integer array holds the number of linear iterations for each Tao iteration
2481: . na    - size of `obj`, `resid`, and `cnorm`
2482: - reset - `PETSC_TRUE` indicates each new minimization resets the history counter to zero,
2483:            else it continues storing new values for new minimizations after the old ones

2485:   Level: intermediate

2487:   Notes:
2488:   If set, `Tao` will fill the given arrays with the indicated
2489:   information at each iteration.  If 'obj','resid','cnorm','lits' are
2490:   *all* `NULL` then space (using size `na`, or 1000 if na is `PETSC_DECIDE` or
2491:   `PETSC_DEFAULT`) is allocated for the history.
2492:   If not all are `NULL`, then only the non-`NULL` information categories
2493:   will be stored, the others will be ignored.

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

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

2501: .seealso: [](ch_tao), `TaoGetConvergenceHistory()`
2502: @*/
2503: PetscErrorCode TaoSetConvergenceHistory(Tao tao, PetscReal obj[], PetscReal resid[], PetscReal cnorm[], PetscInt lits[], PetscInt na, PetscBool reset)
2504: {
2505:   PetscFunctionBegin;
2507:   if (obj) PetscAssertPointer(obj, 2);
2508:   if (resid) PetscAssertPointer(resid, 3);
2509:   if (cnorm) PetscAssertPointer(cnorm, 4);
2510:   if (lits) PetscAssertPointer(lits, 5);

2512:   if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
2513:   if (!obj && !resid && !cnorm && !lits) {
2514:     PetscCall(PetscCalloc4(na, &obj, na, &resid, na, &cnorm, na, &lits));
2515:     tao->hist_malloc = PETSC_TRUE;
2516:   }

2518:   tao->hist_obj   = obj;
2519:   tao->hist_resid = resid;
2520:   tao->hist_cnorm = cnorm;
2521:   tao->hist_lits  = lits;
2522:   tao->hist_max   = na;
2523:   tao->hist_reset = reset;
2524:   tao->hist_len   = 0;
2525:   PetscFunctionReturn(PETSC_SUCCESS);
2526: }

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

2531:   Collective

2533:   Input Parameter:
2534: . tao - the `Tao` context

2536:   Output Parameters:
2537: + obj   - array used to hold objective value history
2538: . resid - array used to hold residual history
2539: . cnorm - array used to hold constraint violation history
2540: . lits  - integer array used to hold linear solver iteration count
2541: - nhist - size of `obj`, `resid`, `cnorm`, and `lits`

2543:   Level: advanced

2545:   Notes:
2546:   This routine must be preceded by calls to `TaoSetConvergenceHistory()`
2547:   and `TaoSolve()`, otherwise it returns useless information.

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

2553:   Fortran Notes:
2554:   The calling sequence is
2555: .vb
2556:    call TaoGetConvergenceHistory(Tao tao, PetscInt nhist, PetscErrorCode ierr)
2557: .ve
2558:   In other words this gets the current number of entries in the history. Access the history through the array you passed to `TaoSetConvergenceHistory()`

2560: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetConvergenceHistory()`
2561: @*/
2562: PetscErrorCode TaoGetConvergenceHistory(Tao tao, PetscReal **obj, PetscReal **resid, PetscReal **cnorm, PetscInt **lits, PetscInt *nhist)
2563: {
2564:   PetscFunctionBegin;
2566:   if (obj) *obj = tao->hist_obj;
2567:   if (cnorm) *cnorm = tao->hist_cnorm;
2568:   if (resid) *resid = tao->hist_resid;
2569:   if (lits) *lits = tao->hist_lits;
2570:   if (nhist) *nhist = tao->hist_len;
2571:   PetscFunctionReturn(PETSC_SUCCESS);
2572: }

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

2577:   Logically Collective

2579:   Input Parameters:
2580: + tao  - the `Tao` context
2581: - usrP - optional user context

2583:   Level: intermediate

2585: .seealso: [](ch_tao), `Tao`, `TaoGetApplicationContext()`
2586: @*/
2587: PetscErrorCode TaoSetApplicationContext(Tao tao, void *usrP)
2588: {
2589:   PetscFunctionBegin;
2591:   tao->user = usrP;
2592:   PetscFunctionReturn(PETSC_SUCCESS);
2593: }

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

2598:   Not Collective

2600:   Input Parameter:
2601: . tao - the `Tao` context

2603:   Output Parameter:
2604: . usrP - user context

2606:   Level: intermediate

2608: .seealso: [](ch_tao), `Tao`, `TaoSetApplicationContext()`
2609: @*/
2610: PetscErrorCode TaoGetApplicationContext(Tao tao, void *usrP)
2611: {
2612:   PetscFunctionBegin;
2614:   PetscAssertPointer(usrP, 2);
2615:   *(void **)usrP = tao->user;
2616:   PetscFunctionReturn(PETSC_SUCCESS);
2617: }

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

2622:   Collective

2624:   Input Parameters:
2625: + tao - the `Tao` context
2626: - M   - matrix that defines the norm

2628:   Level: beginner

2630: .seealso: [](ch_tao), `Tao`, `TaoGetGradientNorm()`, `TaoGradientNorm()`
2631: @*/
2632: PetscErrorCode TaoSetGradientNorm(Tao tao, Mat M)
2633: {
2634:   PetscFunctionBegin;
2637:   PetscCall(PetscObjectReference((PetscObject)M));
2638:   PetscCall(MatDestroy(&tao->gradient_norm));
2639:   PetscCall(VecDestroy(&tao->gradient_norm_tmp));
2640:   tao->gradient_norm = M;
2641:   PetscCall(MatCreateVecs(M, NULL, &tao->gradient_norm_tmp));
2642:   PetscFunctionReturn(PETSC_SUCCESS);
2643: }

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

2648:   Not Collective

2650:   Input Parameter:
2651: . tao - the `Tao` context

2653:   Output Parameter:
2654: . M - gradient norm

2656:   Level: beginner

2658: .seealso: [](ch_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGradientNorm()`
2659: @*/
2660: PetscErrorCode TaoGetGradientNorm(Tao tao, Mat *M)
2661: {
2662:   PetscFunctionBegin;
2664:   PetscAssertPointer(M, 2);
2665:   *M = tao->gradient_norm;
2666:   PetscFunctionReturn(PETSC_SUCCESS);
2667: }

2669: /*@C
2670:   TaoGradientNorm - Compute the norm using the `NormType`, the user has selected

2672:   Collective

2674:   Input Parameters:
2675: + tao      - the `Tao` context
2676: . gradient - the gradient
2677: - type     - the norm type

2679:   Output Parameter:
2680: . gnorm - the gradient norm

2682:   Level: advanced

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

2687:   Developer Notes:
2688:   Should be named `TaoComputeGradientNorm()`.

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

2693: .seealso: [](ch_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGetGradientNorm()`
2694: @*/
2695: PetscErrorCode TaoGradientNorm(Tao tao, Vec gradient, NormType type, PetscReal *gnorm)
2696: {
2697:   PetscFunctionBegin;
2701:   PetscAssertPointer(gnorm, 4);
2702:   if (tao->gradient_norm) {
2703:     PetscScalar gnorms;

2705:     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.");
2706:     PetscCall(MatMult(tao->gradient_norm, gradient, tao->gradient_norm_tmp));
2707:     PetscCall(VecDot(gradient, tao->gradient_norm_tmp, &gnorms));
2708:     *gnorm = PetscRealPart(PetscSqrtScalar(gnorms));
2709:   } else {
2710:     PetscCall(VecNorm(gradient, type, gnorm));
2711:   }
2712:   PetscFunctionReturn(PETSC_SUCCESS);
2713: }

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

2718:   Collective

2720:   Input Parameters:
2721: + comm     - the communicator to share the context
2722: . host     - the name of the X Windows host that will display the monitor
2723: . label    - the label to put at the top of the display window
2724: . x        - the horizontal coordinate of the lower left corner of the window to open
2725: . y        - the vertical coordinate of the lower left corner of the window to open
2726: . m        - the width of the window
2727: . n        - the height of the window
2728: - howoften - how many `Tao` iterations between displaying the monitor information

2730:   Output Parameter:
2731: . ctx - the monitor context

2733:   Options Database Keys:
2734: + -tao_monitor_solution_draw - use `TaoMonitorSolutionDraw()` to monitor the solution
2735: - -tao_draw_solution_initial - show initial guess as well as current solution

2737:   Level: intermediate

2739:   Note:
2740:   The context this creates, along with `TaoMonitorSolutionDraw()`, and `TaoMonitorDrawCtxDestroy()`
2741:   are passed to `TaoMonitorSet()`.

2743: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorDrawCtx()`
2744: @*/
2745: PetscErrorCode TaoMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TaoMonitorDrawCtx *ctx)
2746: {
2747:   PetscFunctionBegin;
2748:   PetscCall(PetscNew(ctx));
2749:   PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
2750:   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
2751:   (*ctx)->howoften = howoften;
2752:   PetscFunctionReturn(PETSC_SUCCESS);
2753: }

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

2758:   Collective

2760:   Input Parameter:
2761: . ictx - the monitor context

2763:   Level: intermediate

2765:   Note:
2766:   This is passed to `TaoMonitorSet()` as the final argument, along with `TaoMonitorSolutionDraw()`, and the context
2767:   obtained with `TaoMonitorDrawCtxCreate()`.

2769: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorSolutionDraw()`
2770: @*/
2771: PetscErrorCode TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx *ictx)
2772: {
2773:   PetscFunctionBegin;
2774:   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
2775:   PetscCall(PetscFree(*ictx));
2776:   PetscFunctionReturn(PETSC_SUCCESS);
2777: }