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 = 0;

  9: PetscLogEvent TAO_Solve;
 10: PetscLogEvent TAO_ResidualEval;
 11: PetscLogEvent TAO_JacobianEval;
 12: PetscLogEvent TAO_ConstraintsEval;

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

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

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

 26:   PetscFunctionBegin;
 27:   if (!snes_ewdummy) PetscFunctionReturn(PETSC_SUCCESS);
 28:   /* populate snes_ewdummy struct values used in KSPPreSolve_SNESEW */
 29:   snes_ewdummy->vec_func = b;
 30:   snes_ewdummy->rtol     = tao->gttol;
 31:   snes_ewdummy->iter     = tao->niter;
 32:   PetscCall(VecNorm(b, NORM_2, &snes_ewdummy->norm));
 33:   PetscCall(KSPPreSolve_SNESEW(ksp, b, x, snes_ewdummy));
 34:   snes_ewdummy->vec_func = NULL;
 35:   PetscFunctionReturn(PETSC_SUCCESS);
 36: }

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

 43:   PetscFunctionBegin;
 44:   if (!snes_ewdummy) PetscFunctionReturn(PETSC_SUCCESS);
 45:   PetscCall(KSPPostSolve_SNESEW(ksp, b, x, snes_ewdummy));
 46:   PetscFunctionReturn(PETSC_SUCCESS);
 47: }

 49: static PetscErrorCode TaoSetUpEW_Private(Tao tao)
 50: {
 51:   SNESKSPEW  *kctx;
 52:   const char *ewprefix;

 54:   PetscFunctionBegin;
 55:   if (!tao->ksp) PetscFunctionReturn(PETSC_SUCCESS);
 56:   if (tao->ksp_ewconv) {
 57:     if (!tao->snes_ewdummy) PetscCall(SNESCreate(PetscObjectComm((PetscObject)tao), &tao->snes_ewdummy));
 58:     tao->snes_ewdummy->ksp_ewconv = PETSC_TRUE;
 59:     PetscCall(KSPSetPreSolve(tao->ksp, KSPPreSolve_TAOEW_Private, tao));
 60:     PetscCall(KSPSetPostSolve(tao->ksp, KSPPostSolve_TAOEW_Private, tao));

 62:     PetscCall(KSPGetOptionsPrefix(tao->ksp, &ewprefix));
 63:     kctx = (SNESKSPEW *)tao->snes_ewdummy->kspconvctx;
 64:     PetscCall(SNESEWSetFromOptions_Private(kctx, PETSC_FALSE, PetscObjectComm((PetscObject)tao), ewprefix));
 65:   } else PetscCall(SNESDestroy(&tao->snes_ewdummy));
 66:   PetscFunctionReturn(PETSC_SUCCESS);
 67: }

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

 73:   Collective

 75:   Input Parameter:
 76: . tao - the `Tao` object

 78:   Level: developer

 80:   Developer Note:
 81:   This is called by all the `TaoCreate_XXX()` routines.

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

101: /*@
102:   TaoCreate - Creates a Tao solver

104:   Collective

106:   Input Parameter:
107: . comm - MPI communicator

109:   Output Parameter:
110: . newtao - the new `Tao` context

112:   Options Database Key:
113: . -tao_type - select which method Tao should use

115:   Level: beginner

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

123:   PetscFunctionBegin;
124:   PetscAssertPointer(newtao, 2);
125:   PetscCall(TaoInitializePackage());
126:   PetscCall(TaoLineSearchInitializePackage());

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

131:   tao->hist_reset = PETSC_TRUE;
132:   tao->term_set   = PETSC_FALSE;

134:   PetscCall(TaoTermCreateCallbacks(tao, &tao->callbacks));
135:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tao->callbacks, "callbacks_"));
136:   PetscCall(TaoTermMappingSetData(&tao->objective_term, NULL, 1.0, tao->callbacks, NULL));
137:   PetscCall(TaoResetStatistics(tao));
138:   *newtao = tao;
139:   PetscFunctionReturn(PETSC_SUCCESS);
140: }

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

145:   Collective

147:   Input Parameter:
148: . tao - the `Tao` context

150:   Level: beginner

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

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

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

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

179:   PetscCall(PetscLogEventBegin(TAO_Solve, tao, 0, 0, 0));
180:   PetscTryTypeMethod(tao, solve);
181:   PetscCall(PetscLogEventEnd(TAO_Solve, tao, 0, 0, 0));

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

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

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

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

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

214:   Collective

216:   Input Parameter:
217: . tao - the `Tao` context

219:   Level: advanced

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

227: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
228: @*/
229: PetscErrorCode TaoSetUp(Tao tao)
230: {
231:   PetscFunctionBegin;
233:   if (tao->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
234:   PetscCall(TaoSetUpEW_Private(tao));
235:   PetscCall(TaoTermMappingSetUp(&tao->objective_term));
236:   if (!tao->solution) PetscCall(TaoTermMappingCreateSolutionVec(&tao->objective_term, &tao->solution));
237:   PetscCheck(tao->solution, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoSetSolution()");
238:   if (tao->uses_gradient && !tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
239:   if (tao->uses_hessian_matrices) {
240:     // TaoSetHessian has been called, but as terms have been added,
241:     // subterms' Hessian and PtAP routines, if needed, have to be created
242:     // TODO Function to set TAOTERMSUM's Hessian.
243:     if (!tao->hessian) {
244:       PetscBool is_defined;

246:       // TAOTERMSUM's Hessian will follow layout and type of first term's Hessian
247:       PetscCall(TaoTermIsCreateHessianMatricesDefined(tao->objective_term.term, &is_defined));
248:       if (is_defined) PetscCall(TaoTermMappingCreateHessianMatrices(&tao->objective_term, &tao->hessian, &tao->hessian_pre));
249:     }
250:     PetscCheck(tao->hessian, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoSetHessian()");
251:   }
252:   PetscTryTypeMethod(tao, setup);
253:   tao->setupcalled = PETSC_TRUE;
254:   PetscFunctionReturn(PETSC_SUCCESS);
255: }

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

260:   Collective

262:   Input Parameter:
263: . tao - the `Tao` context

265:   Level: beginner

267: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
268: @*/
269: PetscErrorCode TaoDestroy(Tao *tao)
270: {
271:   PetscFunctionBegin;
272:   if (!*tao) PetscFunctionReturn(PETSC_SUCCESS);
274:   if (--((PetscObject)*tao)->refct > 0) {
275:     *tao = NULL;
276:     PetscFunctionReturn(PETSC_SUCCESS);
277:   }

279:   PetscTryTypeMethod(*tao, destroy);
280:   PetscCall(TaoTermMappingReset(&(*tao)->objective_term));
281:   PetscCall(VecDestroy(&(*tao)->objective_parameters));
282:   PetscCall(TaoTermDestroy(&(*tao)->callbacks));
283:   PetscCall(KSPDestroy(&(*tao)->ksp));
284:   PetscCall(SNESDestroy(&(*tao)->snes_ewdummy));
285:   PetscCall(TaoLineSearchDestroy(&(*tao)->linesearch));

287:   if ((*tao)->ops->convergencedestroy) {
288:     PetscCall((*(*tao)->ops->convergencedestroy)((*tao)->cnvP));
289:     PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
290:   }
291:   PetscCall(VecDestroy(&(*tao)->solution));
292:   PetscCall(VecDestroy(&(*tao)->gradient));
293:   PetscCall(VecDestroy(&(*tao)->ls_res));

295:   if ((*tao)->gradient_norm) {
296:     PetscCall(PetscObjectDereference((PetscObject)(*tao)->gradient_norm));
297:     PetscCall(VecDestroy(&(*tao)->gradient_norm_tmp));
298:   }

300:   PetscCall(VecDestroy(&(*tao)->XL));
301:   PetscCall(VecDestroy(&(*tao)->XU));
302:   PetscCall(VecDestroy(&(*tao)->IL));
303:   PetscCall(VecDestroy(&(*tao)->IU));
304:   PetscCall(VecDestroy(&(*tao)->DE));
305:   PetscCall(VecDestroy(&(*tao)->DI));
306:   PetscCall(VecDestroy(&(*tao)->constraints));
307:   PetscCall(VecDestroy(&(*tao)->constraints_equality));
308:   PetscCall(VecDestroy(&(*tao)->constraints_inequality));
309:   PetscCall(VecDestroy(&(*tao)->stepdirection));
310:   PetscCall(MatDestroy(&(*tao)->hessian_pre));
311:   PetscCall(MatDestroy(&(*tao)->hessian));
312:   PetscCall(MatDestroy(&(*tao)->ls_jac));
313:   PetscCall(MatDestroy(&(*tao)->ls_jac_pre));
314:   PetscCall(MatDestroy(&(*tao)->jacobian_pre));
315:   PetscCall(MatDestroy(&(*tao)->jacobian));
316:   PetscCall(MatDestroy(&(*tao)->jacobian_state_pre));
317:   PetscCall(MatDestroy(&(*tao)->jacobian_state));
318:   PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
319:   PetscCall(MatDestroy(&(*tao)->jacobian_design));
320:   PetscCall(MatDestroy(&(*tao)->jacobian_equality));
321:   PetscCall(MatDestroy(&(*tao)->jacobian_equality_pre));
322:   PetscCall(MatDestroy(&(*tao)->jacobian_inequality));
323:   PetscCall(MatDestroy(&(*tao)->jacobian_inequality_pre));
324:   PetscCall(ISDestroy(&(*tao)->state_is));
325:   PetscCall(ISDestroy(&(*tao)->design_is));
326:   PetscCall(VecDestroy(&(*tao)->res_weights_v));
327:   PetscCall(TaoMonitorCancel(*tao));
328:   if ((*tao)->hist_malloc) PetscCall(PetscFree4((*tao)->hist_obj, (*tao)->hist_resid, (*tao)->hist_cnorm, (*tao)->hist_lits));
329:   if ((*tao)->res_weights_n) {
330:     PetscCall(PetscFree((*tao)->res_weights_rows));
331:     PetscCall(PetscFree((*tao)->res_weights_cols));
332:     PetscCall(PetscFree((*tao)->res_weights_w));
333:   }
334:   PetscCall(PetscHeaderDestroy(tao));
335:   PetscFunctionReturn(PETSC_SUCCESS);
336: }

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

341:   Logically Collective

343:   Input Parameters:
344: + tao  - Tao context
345: - flag - `PETSC_TRUE` or `PETSC_FALSE`

347:   Level: advanced

349:   Note:
350:   See `SNESKSPSetUseEW()` for customization details.

352: .seealso: [](ch_tao), `Tao`, `SNESKSPSetUseEW()`
353: @*/
354: PetscErrorCode TaoKSPSetUseEW(Tao tao, PetscBool flag)
355: {
356:   PetscFunctionBegin;
359:   tao->ksp_ewconv = flag;
360:   PetscFunctionReturn(PETSC_SUCCESS);
361: }

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

366:   Collective

368:   Input Parameters:
369: + tao     - `Tao` object you wish to monitor
370: . name    - the monitor type one is seeking
371: . help    - message indicating what monitoring is done
372: . manual  - manual page for the monitor
373: - monitor - the monitor function, this must use a `PetscViewerFormat` as its context

375:   Level: developer

377: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
378:           `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
379:           `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
380:           `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
381:           `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
382:           `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
383:           `PetscOptionsFList()`, `PetscOptionsEList()`
384: @*/
385: PetscErrorCode TaoMonitorSetFromOptions(Tao tao, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(Tao, PetscViewerAndFormat *))
386: {
387:   PetscViewer       viewer;
388:   PetscViewerFormat format;
389:   PetscBool         flg;

391:   PetscFunctionBegin;
392:   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)tao), ((PetscObject)tao)->options, ((PetscObject)tao)->prefix, name, &viewer, &format, &flg));
393:   if (flg) {
394:     PetscViewerAndFormat *vf;
395:     char                  interval_key[1024];

397:     PetscCall(PetscSNPrintf(interval_key, sizeof interval_key, "%s_interval", name));
398:     PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
399:     vf->view_interval = 1;
400:     PetscCall(PetscOptionsGetInt(((PetscObject)tao)->options, ((PetscObject)tao)->prefix, interval_key, &vf->view_interval, NULL));

402:     PetscCall(PetscViewerDestroy(&viewer));
403:     PetscCall(TaoMonitorSet(tao, (PetscErrorCode (*)(Tao, PetscCtx))monitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy));
404:   }
405:   PetscFunctionReturn(PETSC_SUCCESS);
406: }

408: /*@
409:   TaoSetFromOptions - Sets various Tao parameters from the options database

411:   Collective

413:   Input Parameter:
414: . tao - the `Tao` solver context

416:   Options Database Keys:
417: + -tao_type type               - The algorithm that Tao uses (lmvm, nls, etc.). See `TAOType`
418: . -tao_gatol gatol             - absolute error tolerance for ||gradient||
419: . -tao_grtol grtol             - relative error tolerance for ||gradient||
420: . -tao_gttol gttol             - reduction of ||gradient|| relative to initial gradient
421: . -tao_max_it max              - sets maximum number of iterations
422: . -tao_max_funcs max           - sets maximum number of function evaluations
423: . -tao_fmin fmin               - stop if function value reaches fmin
424: . -tao_steptol tol             - stop if trust region radius less than `tol`
425: . -tao_trust0 radius           - initial trust region radius
426: . -tao_view_solution           - view the solution at the end of the optimization process
427: . -tao_monitor                 - prints function value and residual norm at each iteration
428: . -tao_monitor_short           - same as `-tao_monitor`, but truncates very small values
429: . -tao_monitor_constraint_norm - prints objective value, gradient, and constraint norm at each iteration
430: . -tao_monitor_globalization   - prints information about the globalization at each iteration
431: . -tao_monitor_solution        - prints solution vector at each iteration
432: . -tao_monitor_ls_residual     - prints least-squares residual vector at each iteration
433: . -tao_monitor_step            - prints step vector at each iteration
434: . -tao_monitor_gradient        - prints gradient vector at each iteration
435: . -tao_monitor_solution_draw   - graphically view solution vector at each iteration
436: . -tao_monitor_step_draw       - graphically view step vector at each iteration
437: . -tao_monitor_gradient_draw   - graphically view gradient at each iteration
438: . -tao_monitor_cancel          - cancels all monitors (except those set with command line)
439: . -tao_fd_gradient             - use gradient computed with finite differences
440: . -tao_fd_hessian              - use hessian computed with finite differences
441: . -tao_mf_hessian              - use matrix-free Hessian computed with finite differences. No `TaoTerm` support
442: . -tao_view                    - prints information about the Tao after solving
443: . -tao_converged_reason        - prints the reason Tao stopped iterating
444: - -tao_add_terms               - takes a comma-separated list of up to 16 options prefixes, a `TaoTerm` will be created for each and added to the objective function

446:   Level: beginner

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

452:   The `-tao_add_terms` option accepts at most 16 prefixes.

454: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
455: @*/
456: PetscErrorCode TaoSetFromOptions(Tao tao)
457: {
458:   TaoType   default_type = TAOLMVM;
459:   char      type[256];
460:   PetscBool flg, found;
461:   MPI_Comm  comm;
462:   PetscReal catol, crtol, gatol, grtol, gttol;

464:   PetscFunctionBegin;
466:   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));

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

470:   PetscObjectOptionsBegin((PetscObject)tao);
471:   /* Check for type from options */
472:   PetscCall(PetscOptionsFList("-tao_type", "Tao Solver type", "TaoSetType", TaoList, default_type, type, 256, &flg));
473:   if (flg) PetscCall(TaoSetType(tao, type));
474:   else if (!((PetscObject)tao)->type_name) PetscCall(TaoSetType(tao, default_type));

476:   /* Tao solvers do not set the prefix, set it here if not yet done
477:      We do it after SetType since solver may have been changed */
478:   if (tao->linesearch) {
479:     const char *prefix;
480:     PetscCall(TaoLineSearchGetOptionsPrefix(tao->linesearch, &prefix));
481:     if (!prefix) PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, ((PetscObject)tao)->prefix));
482:   }

484:   catol = tao->catol;
485:   crtol = tao->crtol;
486:   PetscCall(PetscOptionsReal("-tao_catol", "Stop if constraints violations within", "TaoSetConstraintTolerances", tao->catol, &catol, NULL));
487:   PetscCall(PetscOptionsReal("-tao_crtol", "Stop if relative constraint violations within", "TaoSetConstraintTolerances", tao->crtol, &crtol, NULL));
488:   PetscCall(TaoSetConstraintTolerances(tao, catol, crtol));

490:   gatol = tao->gatol;
491:   grtol = tao->grtol;
492:   gttol = tao->gttol;
493:   PetscCall(PetscOptionsReal("-tao_gatol", "Stop if norm of gradient less than", "TaoSetTolerances", tao->gatol, &gatol, NULL));
494:   PetscCall(PetscOptionsReal("-tao_grtol", "Stop if norm of gradient divided by the function value is less than", "TaoSetTolerances", tao->grtol, &grtol, NULL));
495:   PetscCall(PetscOptionsReal("-tao_gttol", "Stop if the norm of the gradient is less than the norm of the initial gradient times tol", "TaoSetTolerances", tao->gttol, &gttol, NULL));
496:   PetscCall(TaoSetTolerances(tao, gatol, grtol, gttol));

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

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

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

509:   PetscCall(PetscOptionsDeprecated("-tao_solution_monitor", "-tao_monitor_solution", "3.21", NULL));
510:   PetscCall(PetscOptionsDeprecated("-tao_gradient_monitor", "-tao_monitor_gradient", "3.21", NULL));
511:   PetscCall(PetscOptionsDeprecated("-tao_stepdirection_monitor", "-tao_monitor_step", "3.21", NULL));
512:   PetscCall(PetscOptionsDeprecated("-tao_residual_monitor", "-tao_monitor_residual", "3.21", NULL));
513:   PetscCall(PetscOptionsDeprecated("-tao_smonitor", "-tao_monitor_short", "3.21", NULL));
514:   PetscCall(PetscOptionsDeprecated("-tao_cmonitor", "-tao_monitor_constraint_norm", "3.21", NULL));
515:   PetscCall(PetscOptionsDeprecated("-tao_gmonitor", "-tao_monitor_globalization", "3.21", NULL));
516:   PetscCall(PetscOptionsDeprecated("-tao_draw_solution", "-tao_monitor_solution_draw", "3.21", NULL));
517:   PetscCall(PetscOptionsDeprecated("-tao_draw_gradient", "-tao_monitor_gradient_draw", "3.21", NULL));
518:   PetscCall(PetscOptionsDeprecated("-tao_draw_step", "-tao_monitor_step_draw", "3.21", NULL));

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

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

525:   PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_step", "View step vector after each iteration", "TaoMonitorStep", TaoMonitorStep));
526:   PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_residual", "View least-squares residual vector after each iteration", "TaoMonitorResidual", TaoMonitorResidual));
527:   PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor", "Use the default convergence monitor", "TaoMonitorDefault", TaoMonitorDefault));
528:   PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_globalization", "Use the convergence monitor with extra globalization info", "TaoMonitorGlobalization", TaoMonitorGlobalization));
529:   PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_short", "Use the short convergence monitor", "TaoMonitorDefaultShort", TaoMonitorDefaultShort));
530:   PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_constraint_norm", "Use the default convergence monitor with constraint norm", "TaoMonitorConstraintNorm", TaoMonitorConstraintNorm));

532:   flg = PETSC_FALSE;
533:   PetscCall(PetscOptionsDeprecated("-tao_cancelmonitors", "-tao_monitor_cancel", "3.21", NULL));
534:   PetscCall(PetscOptionsBool("-tao_monitor_cancel", "cancel all monitors and call any registered destroy routines", "TaoMonitorCancel", flg, &flg, NULL));
535:   if (flg) PetscCall(TaoMonitorCancel(tao));

537:   flg = PETSC_FALSE;
538:   PetscCall(PetscOptionsBool("-tao_monitor_solution_draw", "Plot solution vector at each iteration", "TaoMonitorSet", flg, &flg, NULL));
539:   if (flg) {
540:     TaoMonitorDrawCtx drawctx;
541:     PetscInt          howoften = 1;
542:     PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &drawctx));
543:     PetscCall(TaoMonitorSet(tao, TaoMonitorSolutionDraw, drawctx, (PetscCtxDestroyFn *)TaoMonitorDrawCtxDestroy));
544:   }

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

550:   flg = PETSC_FALSE;
551:   PetscCall(PetscOptionsBool("-tao_monitor_gradient_draw", "plots gradient at each iteration", "TaoMonitorSet", flg, &flg, NULL));
552:   if (flg) {
553:     TaoMonitorDrawCtx drawctx;
554:     PetscInt          howoften = 1;
555:     PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &drawctx));
556:     PetscCall(TaoMonitorSet(tao, TaoMonitorGradientDraw, drawctx, (PetscCtxDestroyFn *)TaoMonitorDrawCtxDestroy));
557:   }

559:   flg = PETSC_FALSE;
560:   PetscCall(PetscOptionsBool("-tao_fd_gradient", "compute gradient using finite differences", "TaoDefaultComputeGradient", flg, &flg, NULL));
561:   if (flg) PetscCall(TaoTermComputeGradientSetUseFD(tao->objective_term.term, PETSC_TRUE));
562:   flg = PETSC_FALSE;
563:   PetscCall(PetscOptionsBool("-tao_fd_hessian", "compute Hessian using finite differences", "TaoDefaultComputeHessian", flg, &flg, NULL));
564:   if (flg) {
565:     Mat H;

567:     PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &H));
568:     PetscCall(MatSetType(H, MATAIJ));
569:     PetscCall(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE));
570:     PetscCall(MatSetOption(H, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
571:     PetscCall(TaoSetHessian(tao, H, H, TaoDefaultComputeHessian, NULL));
572:     PetscCall(TaoTermComputeHessianSetUseFD(tao->objective_term.term, PETSC_TRUE));
573:     PetscCall(MatDestroy(&H));
574:   }
575:   flg = PETSC_FALSE;
576:   PetscCall(PetscOptionsBool("-tao_mf_hessian", "compute matrix-free Hessian using finite differences", "TaoDefaultComputeHessianMFFD", flg, &flg, NULL));
577:   if (flg) {
578:     PetscBool is_callback;
579:     Mat       H;

581:     // Check that tao has only one TaoTerm with type TAOTERMCALLBACK
582:     PetscCall(PetscObjectTypeCompare((PetscObject)tao->objective_term.term, TAOTERMCALLBACKS, &is_callback));
583:     if (is_callback) {
584:       // Create Hessian via TaoTermCreateHessianMFFD
585:       PetscCall(TaoTermCreateHessianMFFD(tao->objective_term.term, &H));
586:       PetscCall(TaoSetHessian(tao, H, H, TaoDefaultComputeHessianMFFD, NULL));
587:       PetscCall(MatDestroy(&H));
588:     } else {
589:       PetscCall(PetscInfo(tao, "-tao_mf_hessian only works when Tao has a single TAOTERMCALLBACK term. Ignoring.\n"));
590:     }
591:   }
592:   PetscCall(PetscOptionsBool("-tao_recycle_history", "enable recycling/re-using information from the previous TaoSolve() call for some algorithms", "TaoSetRecycleHistory", flg, &flg, &found));
593:   if (found) PetscCall(TaoSetRecycleHistory(tao, flg));
594:   PetscCall(PetscOptionsEnum("-tao_subset_type", "subset type", "", TaoSubSetTypes, (PetscEnum)tao->subset_type, (PetscEnum *)&tao->subset_type, NULL));

596:   if (tao->ksp) {
597:     PetscCall(PetscOptionsBool("-tao_ksp_ew", "Use Eisentat-Walker linear system convergence test", "TaoKSPSetUseEW", tao->ksp_ewconv, &tao->ksp_ewconv, NULL));
598:     PetscCall(TaoKSPSetUseEW(tao, tao->ksp_ewconv));
599:   }

601:   PetscCall(TaoTermSetFromOptions(tao->callbacks));

603:   {
604:     char    *term_prefixes[16];
605:     PetscInt n_terms = PETSC_STATIC_ARRAY_LENGTH(term_prefixes);

607:     PetscCall(PetscOptionsStringArray("-tao_add_terms", "a list of prefixes for terms to add to the Tao objective function", "TaoAddTerm", term_prefixes, &n_terms, NULL));
608:     for (PetscInt i = 0; i < n_terms; i++) {
609:       TaoTerm     term;
610:       const char *prefix;

612:       PetscCall(TaoTermDuplicate(tao->objective_term.term, TAOTERM_DUPLICATE_SIZEONLY, &term));
613:       PetscCall(TaoGetOptionsPrefix(tao, &prefix));
614:       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)term, prefix));
615:       PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)term, term_prefixes[i]));
616:       PetscCall(TaoTermSetFromOptions(term));
617:       PetscCall(TaoAddTerm(tao, term_prefixes[i], 1.0, term, NULL, NULL));
618:       PetscCall(TaoTermDestroy(&term));
619:       PetscCall(PetscFree(term_prefixes[i]));
620:     }
621:   }

623:   if (tao->objective_term.term != tao->callbacks) PetscCall(TaoTermSetFromOptions(tao->objective_term.term));

625:   PetscTryTypeMethod(tao, setfromoptions, PetscOptionsObject);

627:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
628:   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)tao, PetscOptionsObject));
629:   PetscOptionsEnd();

631:   if (tao->linesearch) PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
632:   PetscFunctionReturn(PETSC_SUCCESS);
633: }

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

638:   Collective

640:   Input Parameters:
641: + A    - the  `Tao` context
642: . obj  - Optional object that provides the prefix for the options database
643: - name - command line option

645:   Level: intermediate

647: .seealso: [](ch_tao), `Tao`, `TaoView`, `PetscObjectViewFromOptions()`, `TaoCreate()`
648: @*/
649: PetscErrorCode TaoViewFromOptions(Tao A, PetscObject obj, const char name[])
650: {
651:   PetscFunctionBegin;
653:   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
654:   PetscFunctionReturn(PETSC_SUCCESS);
655: }

657: /*@
658:   TaoView - Prints information about the `Tao` object

660:   Collective

662:   Input Parameters:
663: + tao    - the `Tao` context
664: - viewer - visualization context

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

669:   Level: beginner

671:   Notes:
672:   The available visualization contexts include
673: +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
674: -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
675:   output where only the first processor opens
676:   the file.  All other processors send their
677:   data to the first processor to print.

679:   To view all the `TaoTerm` inside of `Tao`, use `PETSC_VIEWER_ASCII_INFO_DETAIL`,
680:   or pass `-tao_view ::ascii_info_detail` flag

682: .seealso: [](ch_tao), `Tao`, `PetscViewerASCIIOpen()`
683: @*/
684: PetscErrorCode TaoView(Tao tao, PetscViewer viewer)
685: {
686:   PetscBool isascii, isstring;
687:   TaoType   type;

689:   PetscFunctionBegin;
691:   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(((PetscObject)tao)->comm, &viewer));
693:   PetscCheckSameComm(tao, 1, viewer, 2);

695:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
696:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
697:   if (isascii) {
698:     PetscViewerFormat format;

700:     PetscCall(PetscViewerGetFormat(viewer, &format));
701:     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)tao, viewer));

703:     PetscCall(PetscViewerASCIIPushTab(viewer));
704:     PetscTryTypeMethod(tao, view, viewer);
705:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
706:       PetscCall(PetscViewerASCIIPrintf(viewer, "Objective function:\n"));
707:       PetscCall(PetscViewerASCIIPushTab(viewer));
708:       PetscCall(PetscViewerASCIIPrintf(viewer, "Scale (tao_objective_scale): %g\n", (double)tao->objective_term.scale));
709:       PetscCall(PetscViewerASCIIPrintf(viewer, "Function:\n"));
710:       PetscCall(PetscViewerASCIIPushTab(viewer));
711:       PetscCall(TaoTermView(tao->objective_term.term, viewer));
712:       PetscCall(PetscViewerASCIIPopTab(viewer));
713:       if (tao->objective_term.map) {
714:         PetscCall(PetscViewerASCIIPrintf(viewer, "Map:\n"));
715:         PetscCall(PetscViewerASCIIPushTab(viewer));
716:         PetscCall(MatView(tao->objective_term.map, viewer));
717:         PetscCall(PetscViewerASCIIPopTab(viewer));
718:       } else PetscCall(PetscViewerASCIIPrintf(viewer, "Map: unmapped\n"));
719:       PetscCall(PetscViewerASCIIPopTab(viewer));
720:     } else if (tao->num_terms > 0 || tao->term_set) {
721:       if (tao->objective_term.scale == 1.0 && tao->objective_term.map == NULL) {
722:         PetscCall(PetscViewerASCIIPrintf(viewer, "Objective function:\n"));
723:         PetscCall(PetscViewerASCIIPushTab(viewer));
724:         PetscCall(TaoTermView(tao->objective_term.term, viewer));
725:         PetscCall(PetscViewerASCIIPopTab(viewer));
726:       } else {
727:         PetscCall(PetscViewerASCIIPrintf(viewer, "Objective function:\n"));
728:         PetscCall(PetscViewerASCIIPushTab(viewer));
729:         if (tao->objective_term.scale != 1.0) PetscCall(PetscViewerASCIIPrintf(viewer, "Scale: %g\n", (double)tao->objective_term.scale));
730:         PetscCall(PetscViewerASCIIPrintf(viewer, "Function:\n"));
731:         PetscCall(PetscViewerASCIIPushTab(viewer));
732:         PetscCall(TaoTermView(tao->objective_term.term, viewer));
733:         PetscCall(PetscViewerASCIIPopTab(viewer));
734:         if (tao->objective_term.map) {
735:           PetscCall(PetscViewerASCIIPrintf(viewer, "Map:\n"));
736:           PetscCall(PetscViewerASCIIPushTab(viewer));
737:           PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
738:           PetscCall(MatView(tao->objective_term.map, viewer));
739:           PetscCall(PetscViewerPopFormat(viewer));
740:           PetscCall(PetscViewerASCIIPopTab(viewer));
741:         }
742:         PetscCall(PetscViewerASCIIPopTab(viewer));
743:       }
744:     }
745:     if (tao->linesearch) PetscCall(TaoLineSearchView(tao->linesearch, viewer));
746:     if (tao->ksp) {
747:       PetscCall(KSPView(tao->ksp, viewer));
748:       PetscCall(PetscViewerASCIIPrintf(viewer, "total KSP iterations: %" PetscInt_FMT "\n", tao->ksp_tot_its));
749:     }

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

753:     PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: gatol=%g,", (double)tao->gatol));
754:     PetscCall(PetscViewerASCIIPrintf(viewer, " grtol=%g,", (double)tao->grtol));
755:     PetscCall(PetscViewerASCIIPrintf(viewer, " steptol=%g,", (double)tao->steptol));
756:     PetscCall(PetscViewerASCIIPrintf(viewer, " gttol=%g\n", (double)tao->gttol));
757:     PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Function/Gradient:=%g\n", (double)tao->residual));

759:     if (tao->constrained) {
760:       PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances:"));
761:       PetscCall(PetscViewerASCIIPrintf(viewer, " catol=%g,", (double)tao->catol));
762:       PetscCall(PetscViewerASCIIPrintf(viewer, " crtol=%g\n", (double)tao->crtol));
763:       PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Constraints:=%g\n", (double)tao->cnorm));
764:     }

766:     if (tao->trust < tao->steptol) {
767:       PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: steptol=%g\n", (double)tao->steptol));
768:       PetscCall(PetscViewerASCIIPrintf(viewer, "Final trust region radius:=%g\n", (double)tao->trust));
769:     }

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

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

777:     if (tao->objective_term.term->nobj > 0) {
778:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function evaluations=%" PetscInt_FMT ",", tao->objective_term.term->nobj));
779:       if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: unlimited)\n"));
780:       else PetscCall(PetscViewerASCIIPrintf(viewer, "               (max: %" PetscInt_FMT ")\n", tao->max_funcs));
781:     }
782:     if (tao->objective_term.term->ngrad > 0) {
783:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of gradient evaluations=%" PetscInt_FMT ",", tao->objective_term.term->ngrad));
784:       if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: unlimited)\n"));
785:       else PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: %" PetscInt_FMT ")\n", tao->max_funcs));
786:     }
787:     if (tao->objective_term.term->nobjgrad > 0) {
788:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function/gradient evaluations=%" PetscInt_FMT ",", tao->objective_term.term->nobjgrad));
789:       if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, "    (max: unlimited)\n"));
790:       else PetscCall(PetscViewerASCIIPrintf(viewer, "    (max: %" PetscInt_FMT ")\n", tao->max_funcs));
791:     }
792:     if (tao->nres > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of residual evaluations=%" PetscInt_FMT "\n", tao->nres));
793:     if (tao->objective_term.term->nhess > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Hessian evaluations=%" PetscInt_FMT "\n", tao->objective_term.term->nhess));
794:     if (tao->nconstraints > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of constraint function evaluations=%" PetscInt_FMT "\n", tao->nconstraints));
795:     if (tao->njac > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Jacobian evaluations=%" PetscInt_FMT "\n", tao->njac));

797:     if (tao->reason > 0) {
798:       PetscCall(PetscViewerASCIIPrintf(viewer, "Solution converged: "));
799:       switch (tao->reason) {
800:       case TAO_CONVERGED_GATOL:
801:         PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)|| <= gatol\n"));
802:         break;
803:       case TAO_CONVERGED_GRTOL:
804:         PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)||/|f(X)| <= grtol\n"));
805:         break;
806:       case TAO_CONVERGED_GTTOL:
807:         PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)||/||g(X0)|| <= gttol\n"));
808:         break;
809:       case TAO_CONVERGED_STEPTOL:
810:         PetscCall(PetscViewerASCIIPrintf(viewer, " Steptol -- step size small\n"));
811:         break;
812:       case TAO_CONVERGED_MINF:
813:         PetscCall(PetscViewerASCIIPrintf(viewer, " Minf --  f < fmin\n"));
814:         break;
815:       case TAO_CONVERGED_USER:
816:         PetscCall(PetscViewerASCIIPrintf(viewer, " User Terminated\n"));
817:         break;
818:       default:
819:         PetscCall(PetscViewerASCIIPrintf(viewer, " %d\n", tao->reason));
820:         break;
821:       }
822:     } else if (tao->reason == TAO_CONTINUE_ITERATING) {
823:       PetscCall(PetscViewerASCIIPrintf(viewer, "Solver never run\n"));
824:     } else {
825:       PetscCall(PetscViewerASCIIPrintf(viewer, "Solver failed: "));
826:       switch (tao->reason) {
827:       case TAO_DIVERGED_MAXITS:
828:         PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum Iterations\n"));
829:         break;
830:       case TAO_DIVERGED_NAN:
831:         PetscCall(PetscViewerASCIIPrintf(viewer, " NAN or infinity encountered\n"));
832:         break;
833:       case TAO_DIVERGED_MAXFCN:
834:         PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum Function Evaluations\n"));
835:         break;
836:       case TAO_DIVERGED_LS_FAILURE:
837:         PetscCall(PetscViewerASCIIPrintf(viewer, " Line Search Failure\n"));
838:         break;
839:       case TAO_DIVERGED_TR_REDUCTION:
840:         PetscCall(PetscViewerASCIIPrintf(viewer, " Trust Region too small\n"));
841:         break;
842:       case TAO_DIVERGED_USER:
843:         PetscCall(PetscViewerASCIIPrintf(viewer, " User Terminated\n"));
844:         break;
845:       default:
846:         PetscCall(PetscViewerASCIIPrintf(viewer, " %d\n", tao->reason));
847:         break;
848:       }
849:     }
850:     PetscCall(PetscViewerASCIIPopTab(viewer));
851:   } else if (isstring) {
852:     PetscCall(TaoGetType(tao, &type));
853:     PetscCall(PetscViewerStringSPrintf(viewer, " %-3.3s", type));
854:   }
855:   PetscFunctionReturn(PETSC_SUCCESS);
856: }

858: /*@
859:   TaoSetRecycleHistory - Sets the boolean flag to enable/disable re-using
860:   iterate information from the previous `TaoSolve()`. This feature is disabled by
861:   default.

863:   Logically Collective

865:   Input Parameters:
866: + tao     - the `Tao` context
867: - recycle - boolean flag

869:   Options Database Key:
870: . -tao_recycle_history (true|false) - reuse the history

872:   Level: intermediate

874:   Notes:
875:   For conjugate gradient methods (`TAOBNCG`), this re-uses the latest search direction
876:   from the previous `TaoSolve()` call when computing the first search direction in a
877:   new solution. By default, CG methods set the first search direction to the
878:   negative gradient.

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

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

887: .seealso: [](ch_tao), `Tao`, `TaoGetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
888: @*/
889: PetscErrorCode TaoSetRecycleHistory(Tao tao, PetscBool recycle)
890: {
891:   PetscFunctionBegin;
894:   tao->recycle = recycle;
895:   PetscFunctionReturn(PETSC_SUCCESS);
896: }

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

902:   Logically Collective

904:   Input Parameter:
905: . tao - the `Tao` context

907:   Output Parameter:
908: . recycle - boolean flag

910:   Level: intermediate

912: .seealso: [](ch_tao), `Tao`, `TaoSetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
913: @*/
914: PetscErrorCode TaoGetRecycleHistory(Tao tao, PetscBool *recycle)
915: {
916:   PetscFunctionBegin;
918:   PetscAssertPointer(recycle, 2);
919:   *recycle = tao->recycle;
920:   PetscFunctionReturn(PETSC_SUCCESS);
921: }

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

926:   Logically Collective

928:   Input Parameters:
929: + tao   - the `Tao` context
930: . gatol - stop if norm of gradient is less than this
931: . grtol - stop if relative norm of gradient is less than this
932: - gttol - stop if norm of gradient is reduced by this factor

934:   Options Database Keys:
935: + -tao_gatol gatol - Sets gatol
936: . -tao_grtol grtol - Sets grtol
937: - -tao_gttol gttol - Sets gttol

939:   Stopping Criteria\:
940: .vb
941:   ||g(X)||                            <= gatol
942:   ||g(X)|| / |f(X)|                   <= grtol
943:   ||g(X)|| / ||g(X0)||                <= gttol
944: .ve

946:   Level: beginner

948:   Notes:
949:   Use `PETSC_CURRENT` to leave one or more tolerances unchanged.

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

953:   Fortran Note:
954:   Use `PETSC_CURRENT_REAL` or `PETSC_DETERMINE_REAL`

956: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`
957: @*/
958: PetscErrorCode TaoSetTolerances(Tao tao, PetscReal gatol, PetscReal grtol, PetscReal gttol)
959: {
960:   PetscFunctionBegin;

966:   if (gatol == (PetscReal)PETSC_DETERMINE) {
967:     tao->gatol = tao->default_gatol;
968:   } else if (gatol != (PetscReal)PETSC_CURRENT) {
969:     PetscCheck(gatol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative gatol not allowed");
970:     tao->gatol = gatol;
971:   }

973:   if (grtol == (PetscReal)PETSC_DETERMINE) {
974:     tao->grtol = tao->default_grtol;
975:   } else if (grtol != (PetscReal)PETSC_CURRENT) {
976:     PetscCheck(grtol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative grtol not allowed");
977:     tao->grtol = grtol;
978:   }

980:   if (gttol == (PetscReal)PETSC_DETERMINE) {
981:     tao->gttol = tao->default_gttol;
982:   } else if (gttol != (PetscReal)PETSC_CURRENT) {
983:     PetscCheck(gttol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative gttol not allowed");
984:     tao->gttol = gttol;
985:   }
986:   PetscFunctionReturn(PETSC_SUCCESS);
987: }

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

992:   Logically Collective

994:   Input Parameters:
995: + tao   - the `Tao` context
996: . catol - absolute constraint tolerance, constraint norm must be less than `catol` for used for `gatol` convergence criteria
997: - crtol - relative constraint tolerance, constraint norm must be less than `crtol` for used for `gatol`, `gttol` convergence criteria

999:   Options Database Keys:
1000: + -tao_catol catol - Sets catol
1001: - -tao_crtol crtol - Sets crtol

1003:   Level: intermediate

1005:   Notes:
1006:   Use `PETSC_CURRENT` to leave one or tolerance unchanged.

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

1010:   Fortran Note:
1011:   Use `PETSC_CURRENT_REAL` or `PETSC_DETERMINE_REAL`

1013: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`, `TaoGetConstraintTolerances()`, `TaoSetTolerances()`
1014: @*/
1015: PetscErrorCode TaoSetConstraintTolerances(Tao tao, PetscReal catol, PetscReal crtol)
1016: {
1017:   PetscFunctionBegin;

1022:   if (catol == (PetscReal)PETSC_DETERMINE) {
1023:     tao->catol = tao->default_catol;
1024:   } else if (catol != (PetscReal)PETSC_CURRENT) {
1025:     PetscCheck(catol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative catol not allowed");
1026:     tao->catol = catol;
1027:   }

1029:   if (crtol == (PetscReal)PETSC_DETERMINE) {
1030:     tao->crtol = tao->default_crtol;
1031:   } else if (crtol != (PetscReal)PETSC_CURRENT) {
1032:     PetscCheck(crtol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative crtol not allowed");
1033:     tao->crtol = crtol;
1034:   }
1035:   PetscFunctionReturn(PETSC_SUCCESS);
1036: }

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

1041:   Not Collective

1043:   Input Parameter:
1044: . tao - the `Tao` context

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

1050:   Level: intermediate

1052: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`, `TaoSetTolerances()`, `TaoSetConstraintTolerances()`
1053: @*/
1054: PetscErrorCode TaoGetConstraintTolerances(Tao tao, PetscReal *catol, PetscReal *crtol)
1055: {
1056:   PetscFunctionBegin;
1058:   if (catol) *catol = tao->catol;
1059:   if (crtol) *crtol = tao->crtol;
1060:   PetscFunctionReturn(PETSC_SUCCESS);
1061: }

1063: /*@
1064:   TaoSetFunctionLowerBound - Sets a bound on the solution objective value.
1065:   When an approximate solution with an objective value below this number
1066:   has been found, the solver will terminate.

1068:   Logically Collective

1070:   Input Parameters:
1071: + tao  - the Tao solver context
1072: - fmin - the tolerance

1074:   Options Database Key:
1075: . -tao_fmin fmin - sets the minimum function value

1077:   Level: intermediate

1079: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetTolerances()`
1080: @*/
1081: PetscErrorCode TaoSetFunctionLowerBound(Tao tao, PetscReal fmin)
1082: {
1083:   PetscFunctionBegin;
1086:   tao->fmin = fmin;
1087:   PetscFunctionReturn(PETSC_SUCCESS);
1088: }

1090: /*@
1091:   TaoGetFunctionLowerBound - Gets the bound on the solution objective value.
1092:   When an approximate solution with an objective value below this number
1093:   has been found, the solver will terminate.

1095:   Not Collective

1097:   Input Parameter:
1098: . tao - the `Tao` solver context

1100:   Output Parameter:
1101: . fmin - the minimum function value

1103:   Level: intermediate

1105: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetFunctionLowerBound()`
1106: @*/
1107: PetscErrorCode TaoGetFunctionLowerBound(Tao tao, PetscReal *fmin)
1108: {
1109:   PetscFunctionBegin;
1111:   PetscAssertPointer(fmin, 2);
1112:   *fmin = tao->fmin;
1113:   PetscFunctionReturn(PETSC_SUCCESS);
1114: }

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

1119:   Logically Collective

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

1125:   Options Database Key:
1126: . -tao_max_funcs nfcn - sets the maximum number of function evaluations

1128:   Level: intermediate

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

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

1136: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoSetMaximumIterations()`
1137: @*/
1138: PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao tao, PetscInt nfcn)
1139: {
1140:   PetscFunctionBegin;
1143:   if (nfcn == PETSC_DETERMINE) {
1144:     tao->max_funcs = tao->default_max_funcs;
1145:   } else if (nfcn == PETSC_UNLIMITED || nfcn < 0) {
1146:     tao->max_funcs = PETSC_UNLIMITED;
1147:   } else {
1148:     PetscCheck(nfcn >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of function evaluations  must be positive");
1149:     tao->max_funcs = nfcn;
1150:   }
1151:   PetscFunctionReturn(PETSC_SUCCESS);
1152: }

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

1157:   Logically Collective

1159:   Input Parameter:
1160: . tao - the `Tao` solver context

1162:   Output Parameter:
1163: . nfcn - the maximum number of function evaluations

1165:   Level: intermediate

1167: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1168: @*/
1169: PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao tao, PetscInt *nfcn)
1170: {
1171:   PetscFunctionBegin;
1173:   PetscAssertPointer(nfcn, 2);
1174:   *nfcn = tao->max_funcs;
1175:   PetscFunctionReturn(PETSC_SUCCESS);
1176: }

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

1181:   Not Collective

1183:   Input Parameter:
1184: . tao - the `Tao` solver context

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

1189:   Level: intermediate

1191: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1192: @*/
1193: PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao tao, PetscInt *nfuncs)
1194: {
1195:   PetscFunctionBegin;
1197:   PetscAssertPointer(nfuncs, 2);
1198:   *nfuncs = PetscMax(tao->objective_term.term->nobj, tao->objective_term.term->nobjgrad);
1199:   PetscFunctionReturn(PETSC_SUCCESS);
1200: }

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

1205:   Logically Collective

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

1211:   Options Database Key:
1212: . -tao_max_it its - sets the maximum number of iterations

1214:   Level: intermediate

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

1219:   Developer Note:
1220:   Also accepts the deprecated negative values to indicate no limit

1222: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoSetMaximumFunctionEvaluations()`
1223: @*/
1224: PetscErrorCode TaoSetMaximumIterations(Tao tao, PetscInt maxits)
1225: {
1226:   PetscFunctionBegin;
1229:   if (maxits == PETSC_DETERMINE) {
1230:     tao->max_it = tao->default_max_it;
1231:   } else if (maxits == PETSC_UNLIMITED) {
1232:     tao->max_it = PETSC_INT_MAX;
1233:   } else {
1234:     PetscCheck(maxits > 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations must be positive");
1235:     tao->max_it = maxits;
1236:   }
1237:   PetscFunctionReturn(PETSC_SUCCESS);
1238: }

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

1243:   Not Collective

1245:   Input Parameter:
1246: . tao - the `Tao` solver context

1248:   Output Parameter:
1249: . maxits - the maximum number of iterates

1251:   Level: intermediate

1253: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumIterations()`, `TaoGetMaximumFunctionEvaluations()`
1254: @*/
1255: PetscErrorCode TaoGetMaximumIterations(Tao tao, PetscInt *maxits)
1256: {
1257:   PetscFunctionBegin;
1259:   PetscAssertPointer(maxits, 2);
1260:   *maxits = tao->max_it;
1261:   PetscFunctionReturn(PETSC_SUCCESS);
1262: }

1264: /*@
1265:   TaoSetInitialTrustRegionRadius - Sets the initial trust region radius.

1267:   Logically Collective

1269:   Input Parameters:
1270: + tao    - a `Tao` optimization solver
1271: - radius - the trust region radius

1273:   Options Database Key:
1274: . -tao_trust0 radius - sets initial trust region radius

1276:   Level: intermediate

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

1281: .seealso: [](ch_tao), `Tao`, `TaoGetTrustRegionRadius()`, `TaoSetTrustRegionTolerance()`, `TAONTR`
1282: @*/
1283: PetscErrorCode TaoSetInitialTrustRegionRadius(Tao tao, PetscReal radius)
1284: {
1285:   PetscFunctionBegin;
1288:   if (radius == PETSC_DETERMINE) {
1289:     tao->trust0 = tao->default_trust0;
1290:   } else {
1291:     PetscCheck(radius > 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Radius must be positive");
1292:     tao->trust0 = radius;
1293:   }
1294:   PetscFunctionReturn(PETSC_SUCCESS);
1295: }

1297: /*@
1298:   TaoGetInitialTrustRegionRadius - Gets the initial trust region radius.

1300:   Not Collective

1302:   Input Parameter:
1303: . tao - a `Tao` optimization solver

1305:   Output Parameter:
1306: . radius - the trust region radius

1308:   Level: intermediate

1310: .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetCurrentTrustRegionRadius()`, `TAONTR`
1311: @*/
1312: PetscErrorCode TaoGetInitialTrustRegionRadius(Tao tao, PetscReal *radius)
1313: {
1314:   PetscFunctionBegin;
1316:   PetscAssertPointer(radius, 2);
1317:   *radius = tao->trust0;
1318:   PetscFunctionReturn(PETSC_SUCCESS);
1319: }

1321: /*@
1322:   TaoGetCurrentTrustRegionRadius - Gets the current trust region radius.

1324:   Not Collective

1326:   Input Parameter:
1327: . tao - a `Tao` optimization solver

1329:   Output Parameter:
1330: . radius - the trust region radius

1332:   Level: intermediate

1334: .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetInitialTrustRegionRadius()`, `TAONTR`
1335: @*/
1336: PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao tao, PetscReal *radius)
1337: {
1338:   PetscFunctionBegin;
1340:   PetscAssertPointer(radius, 2);
1341:   *radius = tao->trust;
1342:   PetscFunctionReturn(PETSC_SUCCESS);
1343: }

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

1348:   Not Collective

1350:   Input Parameter:
1351: . tao - the `Tao` context

1353:   Output Parameters:
1354: + gatol - stop if norm of gradient is less than this
1355: . grtol - stop if relative norm of gradient is less than this
1356: - gttol - stop if norm of gradient is reduced by a this factor

1358:   Level: intermediate

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

1363: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`
1364: @*/
1365: PetscErrorCode TaoGetTolerances(Tao tao, PetscReal *gatol, PetscReal *grtol, PetscReal *gttol)
1366: {
1367:   PetscFunctionBegin;
1369:   if (gatol) *gatol = tao->gatol;
1370:   if (grtol) *grtol = tao->grtol;
1371:   if (gttol) *gttol = tao->gttol;
1372:   PetscFunctionReturn(PETSC_SUCCESS);
1373: }

1375: /*@
1376:   TaoGetKSP - Gets the linear solver used by the optimization solver.

1378:   Not Collective

1380:   Input Parameter:
1381: . tao - the `Tao` solver

1383:   Output Parameter:
1384: . ksp - the `KSP` linear solver used in the optimization solver

1386:   Level: intermediate

1388: .seealso: [](ch_tao), `Tao`, `KSP`
1389: @*/
1390: PetscErrorCode TaoGetKSP(Tao tao, KSP *ksp)
1391: {
1392:   PetscFunctionBegin;
1394:   PetscAssertPointer(ksp, 2);
1395:   *ksp = tao->ksp;
1396:   PetscFunctionReturn(PETSC_SUCCESS);
1397: }

1399: /*@
1400:   TaoGetLinearSolveIterations - Gets the total number of linear iterations
1401:   used by the `Tao` solver

1403:   Not Collective

1405:   Input Parameter:
1406: . tao - the `Tao` context

1408:   Output Parameter:
1409: . lits - number of linear iterations

1411:   Level: intermediate

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

1416: .seealso: [](ch_tao), `Tao`, `TaoGetKSP()`
1417: @*/
1418: PetscErrorCode TaoGetLinearSolveIterations(Tao tao, PetscInt *lits)
1419: {
1420:   PetscFunctionBegin;
1422:   PetscAssertPointer(lits, 2);
1423:   *lits = tao->ksp_tot_its;
1424:   PetscFunctionReturn(PETSC_SUCCESS);
1425: }

1427: /*@
1428:   TaoGetLineSearch - Gets the line search used by the optimization solver.

1430:   Not Collective

1432:   Input Parameter:
1433: . tao - the `Tao` solver

1435:   Output Parameter:
1436: . ls - the line search used in the optimization solver

1438:   Level: intermediate

1440: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`
1441: @*/
1442: PetscErrorCode TaoGetLineSearch(Tao tao, TaoLineSearch *ls)
1443: {
1444:   PetscFunctionBegin;
1446:   PetscAssertPointer(ls, 2);
1447:   *ls = tao->linesearch;
1448:   PetscFunctionReturn(PETSC_SUCCESS);
1449: }

1451: /*@
1452:   TaoAddLineSearchCounts - Adds the number of function evaluations spent
1453:   in the line search to the running total.

1455:   Input Parameters:
1456: . tao - the `Tao` solver

1458:   Level: developer

1460: .seealso: [](ch_tao), `Tao`, `TaoGetLineSearch()`, `TaoLineSearchApply()`
1461: @*/
1462: PetscErrorCode TaoAddLineSearchCounts(Tao tao)
1463: {
1464:   PetscBool flg;
1465:   PetscInt  nfeval, ngeval, nfgeval;

1467:   PetscFunctionBegin;
1469:   if (tao->linesearch) {
1470:     PetscCall(TaoLineSearchIsUsingTaoRoutines(tao->linesearch, &flg));
1471:     if (!flg) {
1472:       PetscCall(TaoLineSearchGetNumberFunctionEvaluations(tao->linesearch, &nfeval, &ngeval, &nfgeval));
1473:       tao->objective_term.term->nobj += nfeval;
1474:       tao->objective_term.term->ngrad += ngeval;
1475:       tao->objective_term.term->nobjgrad += nfgeval;
1476:     }
1477:   }
1478:   PetscFunctionReturn(PETSC_SUCCESS);
1479: }

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

1484:   Not Collective

1486:   Input Parameter:
1487: . tao - the `Tao` context

1489:   Output Parameter:
1490: . X - the current solution

1492:   Level: intermediate

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

1497: .seealso: [](ch_tao), `Tao`, `TaoSetSolution()`, `TaoSolve()`
1498: @*/
1499: PetscErrorCode TaoGetSolution(Tao tao, Vec *X)
1500: {
1501:   PetscFunctionBegin;
1503:   PetscAssertPointer(X, 2);
1504:   *X = tao->solution;
1505:   PetscFunctionReturn(PETSC_SUCCESS);
1506: }

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

1513:   Collective

1515:   Input Parameter:
1516: . tao - the `Tao` context

1518:   Level: developer

1520:   Note:
1521:   This function does not reset the statistics of internal `TaoTerm`

1523: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
1524: @*/
1525: PetscErrorCode TaoResetStatistics(Tao tao)
1526: {
1527:   PetscFunctionBegin;
1529:   tao->niter        = 0;
1530:   tao->nres         = 0;
1531:   tao->njac         = 0;
1532:   tao->nconstraints = 0;
1533:   tao->ksp_its      = 0;
1534:   tao->ksp_tot_its  = 0;
1535:   tao->reason       = TAO_CONTINUE_ITERATING;
1536:   tao->residual     = 0.0;
1537:   tao->cnorm        = 0.0;
1538:   tao->step         = 0.0;
1539:   tao->lsflag       = PETSC_FALSE;
1540:   if (tao->hist_reset) tao->hist_len = 0;
1541:   PetscFunctionReturn(PETSC_SUCCESS);
1542: }

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

1549:   Logically Collective

1551:   Input Parameters:
1552: + tao  - The `Tao` solver
1553: . func - The function
1554: - ctx  - The update function context

1556:   Calling sequence of `func`:
1557: + tao - The optimizer context
1558: . it  - The current iteration index
1559: - ctx - The update context

1561:   Level: advanced

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

1567: .seealso: [](ch_tao), `Tao`, `TaoSolve()`
1568: @*/
1569: PetscErrorCode TaoSetUpdate(Tao tao, PetscErrorCode (*func)(Tao tao, PetscInt it, PetscCtx ctx), PetscCtx ctx)
1570: {
1571:   PetscFunctionBegin;
1573:   tao->ops->update = func;
1574:   tao->user_update = ctx;
1575:   PetscFunctionReturn(PETSC_SUCCESS);
1576: }

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

1583:   Logically Collective

1585:   Input Parameters:
1586: + tao  - the `Tao` object
1587: . conv - the routine to test for convergence
1588: - ctx  - [optional] context for private data for the convergence routine (may be `NULL`)

1590:   Calling sequence of `conv`:
1591: + tao - the `Tao` object
1592: - ctx - [optional] convergence context

1594:   Level: advanced

1596:   Note:
1597:   The new convergence testing routine should call `TaoSetConvergedReason()`.

1599: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetConvergedReason()`, `TaoGetSolutionStatus()`, `TaoGetTolerances()`, `TaoMonitorSet()`
1600: @*/
1601: PetscErrorCode TaoSetConvergenceTest(Tao tao, PetscErrorCode (*conv)(Tao tao, PetscCtx ctx), PetscCtx ctx)
1602: {
1603:   PetscFunctionBegin;
1605:   tao->ops->convergencetest = conv;
1606:   tao->cnvP                 = ctx;
1607:   PetscFunctionReturn(PETSC_SUCCESS);
1608: }

1610: /*@C
1611:   TaoMonitorSet - Sets an additional function that is to be used at every
1612:   iteration of the solver to display the iteration's
1613:   progress.

1615:   Logically Collective

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

1623:   Calling sequence of `func`:
1624: + tao - the `Tao` solver context
1625: - ctx - [optional] monitoring context

1627:   Level: intermediate

1629:   Notes:
1630:   See `TaoSetFromOptions()` for a monitoring options.

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

1636:   Fortran Notes:
1637:   Only one monitor function may be set

1639: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoMonitorDefault()`, `TaoMonitorCancel()`, `TaoView()`, `PetscCtxDestroyFn`
1640: @*/
1641: PetscErrorCode TaoMonitorSet(Tao tao, PetscErrorCode (*func)(Tao tao, PetscCtx ctx), PetscCtx ctx, PetscCtxDestroyFn *dest)
1642: {
1643:   PetscFunctionBegin;
1645:   PetscCheck(tao->numbermonitors < MAXTAOMONITORS, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Cannot attach another monitor -- max=%d", MAXTAOMONITORS);
1646:   for (PetscInt i = 0; i < tao->numbermonitors; i++) {
1647:     PetscBool identical;

1649:     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))(PetscVoidFn *)func, ctx, dest, (PetscErrorCode (*)(void))(PetscVoidFn *)tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i], &identical));
1650:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
1651:   }
1652:   tao->monitor[tao->numbermonitors]        = func;
1653:   tao->monitorcontext[tao->numbermonitors] = ctx;
1654:   tao->monitordestroy[tao->numbermonitors] = dest;
1655:   ++tao->numbermonitors;
1656:   PetscFunctionReturn(PETSC_SUCCESS);
1657: }

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

1662:   Logically Collective

1664:   Input Parameter:
1665: . tao - the `Tao` solver context

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

1672:   Level: advanced

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

1677: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1678: @*/
1679: PetscErrorCode TaoMonitorCancel(Tao tao)
1680: {
1681:   PetscInt i;

1683:   PetscFunctionBegin;
1685:   for (i = 0; i < tao->numbermonitors; i++) {
1686:     if (tao->monitordestroy[i]) PetscCall((*tao->monitordestroy[i])(&tao->monitorcontext[i]));
1687:   }
1688:   tao->numbermonitors = 0;
1689:   PetscFunctionReturn(PETSC_SUCCESS);
1690: }

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

1695:   Collective

1697:   Input Parameters:
1698: + tao - the `Tao` context
1699: - vf  - `PetscViewerAndFormat` context

1701:   Options Database Key:
1702: . -tao_monitor - turn on default monitoring

1704:   Level: advanced

1706:   Note:
1707:   This monitor prints the function value and gradient
1708:   norm at each iteration.

1710: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1711: @*/
1712: PetscErrorCode TaoMonitorDefault(Tao tao, PetscViewerAndFormat *vf)
1713: {
1714:   PetscViewer viewer = vf->viewer;
1715:   PetscBool   isascii;
1716:   PetscInt    tabs;

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

1722:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1723:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1724:   if (isascii) {
1725:     PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));

1727:     PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1728:     if (tao->niter == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
1729:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Iteration information for %s solve.\n", ((PetscObject)tao)->prefix));
1730:       tao->header_printed = PETSC_TRUE;
1731:     }
1732:     PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " TAO,", tao->niter));
1733:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Function value: %g,", (double)tao->fc));
1734:     if (tao->residual >= PETSC_INFINITY) {
1735:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: infinity \n"));
1736:     } else {
1737:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g \n", (double)tao->residual));
1738:     }
1739:     PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1740:   }
1741:   PetscCall(PetscViewerPopFormat(viewer));
1742:   PetscFunctionReturn(PETSC_SUCCESS);
1743: }

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

1748:   Collective

1750:   Input Parameters:
1751: + tao - the `Tao` context
1752: - vf  - `PetscViewerAndFormat` context

1754:   Options Database Key:
1755: . -tao_monitor_globalization - turn on monitoring with globalization information

1757:   Level: advanced

1759:   Note:
1760:   This monitor prints the function value and gradient norm at each
1761:   iteration, as well as the step size and trust radius. Note that the
1762:   step size and trust radius may be the same for some algorithms.

1764: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1765: @*/
1766: PetscErrorCode TaoMonitorGlobalization(Tao tao, PetscViewerAndFormat *vf)
1767: {
1768:   PetscViewer viewer = vf->viewer;
1769:   PetscBool   isascii;
1770:   PetscInt    tabs;

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

1776:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1777:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1778:   if (isascii) {
1779:     PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1780:     PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1781:     if (tao->niter == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
1782:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Iteration information for %s solve.\n", ((PetscObject)tao)->prefix));
1783:       tao->header_printed = PETSC_TRUE;
1784:     }
1785:     PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " TAO,", tao->niter));
1786:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Function value: %g,", (double)tao->fc));
1787:     if (tao->residual >= PETSC_INFINITY) {
1788:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: Inf,"));
1789:     } else {
1790:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g,", (double)tao->residual));
1791:     }
1792:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Step: %g,  Trust: %g\n", (double)tao->step, (double)tao->trust));
1793:     PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1794:   }
1795:   PetscCall(PetscViewerPopFormat(viewer));
1796:   PetscFunctionReturn(PETSC_SUCCESS);
1797: }

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

1802:   Collective

1804:   Input Parameters:
1805: + tao - the `Tao` context
1806: - vf  - `PetscViewerAndFormat` context

1808:   Options Database Key:
1809: . -tao_monitor_short - turn on default short monitoring

1811:   Level: advanced

1813:   Note:
1814:   Same as `TaoMonitorDefault()` except
1815:   it prints fewer digits of the residual as the residual gets smaller.
1816:   This is because the later digits are meaningless and are often
1817:   different on different machines; by using this routine different
1818:   machines will usually generate the same output.

1820: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1821: @*/
1822: PetscErrorCode TaoMonitorDefaultShort(Tao tao, PetscViewerAndFormat *vf)
1823: {
1824:   PetscViewer viewer = vf->viewer;
1825:   PetscBool   isascii;
1826:   PetscInt    tabs;
1827:   PetscReal   gnorm;

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

1833:   gnorm = tao->residual;
1834:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1835:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1836:   if (isascii) {
1837:     PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1838:     PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1839:     PetscCall(PetscViewerASCIIPrintf(viewer, "iter = %3" PetscInt_FMT ",", tao->niter));
1840:     PetscCall(PetscViewerASCIIPrintf(viewer, " Function value %g,", (double)tao->fc));
1841:     if (gnorm >= PETSC_INFINITY) {
1842:       PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: infinity \n"));
1843:     } else if (gnorm > 1.e-6) {
1844:       PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: %g \n", (double)gnorm));
1845:     } else if (gnorm > 1.e-11) {
1846:       PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: < 1.0e-6 \n"));
1847:     } else {
1848:       PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: < 1.0e-11 \n"));
1849:     }
1850:     PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1851:   }
1852:   PetscCall(PetscViewerPopFormat(viewer));
1853:   PetscFunctionReturn(PETSC_SUCCESS);
1854: }

1856: /*@
1857:   TaoMonitorConstraintNorm - same as `TaoMonitorDefault()` except
1858:   it prints the norm of the constraint function.

1860:   Collective

1862:   Input Parameters:
1863: + tao - the `Tao` context
1864: - vf  - `PetscViewerAndFormat` context

1866:   Options Database Key:
1867: . -tao_monitor_constraint_norm - monitor the constraints

1869:   Level: advanced

1871: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1872: @*/
1873: PetscErrorCode TaoMonitorConstraintNorm(Tao tao, PetscViewerAndFormat *vf)
1874: {
1875:   PetscViewer viewer = vf->viewer;
1876:   PetscBool   isascii;
1877:   PetscInt    tabs;

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

1883:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1884:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1885:   if (isascii) {
1886:     PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1887:     PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1888:     PetscCall(PetscViewerASCIIPrintf(viewer, "iter = %" PetscInt_FMT ",", tao->niter));
1889:     PetscCall(PetscViewerASCIIPrintf(viewer, " Function value: %g,", (double)tao->fc));
1890:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g ", (double)tao->residual));
1891:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Constraint: %g \n", (double)tao->cnorm));
1892:     PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1893:   }
1894:   PetscCall(PetscViewerPopFormat(viewer));
1895:   PetscFunctionReturn(PETSC_SUCCESS);
1896: }

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

1901:   Collective

1903:   Input Parameters:
1904: + tao - the `Tao` context
1905: - vf  - `PetscViewerAndFormat` context

1907:   Options Database Key:
1908: . -tao_monitor_solution - view the solution

1910:   Level: advanced

1912: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1913: @*/
1914: PetscErrorCode TaoMonitorSolution(Tao tao, PetscViewerAndFormat *vf)
1915: {
1916:   PetscFunctionBegin;
1918:   if (vf->view_interval > 0 && tao->niter % vf->view_interval) PetscFunctionReturn(PETSC_SUCCESS);
1919:   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
1920:   PetscCall(VecView(tao->solution, vf->viewer));
1921:   PetscCall(PetscViewerPopFormat(vf->viewer));
1922:   PetscFunctionReturn(PETSC_SUCCESS);
1923: }

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

1928:   Collective

1930:   Input Parameters:
1931: + tao - the `Tao` context
1932: - vf  - `PetscViewerAndFormat` context

1934:   Options Database Key:
1935: . -tao_monitor_gradient - view the gradient at each iteration

1937:   Level: advanced

1939: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1940: @*/
1941: PetscErrorCode TaoMonitorGradient(Tao tao, PetscViewerAndFormat *vf)
1942: {
1943:   PetscFunctionBegin;
1945:   if (vf->view_interval > 0 && tao->niter % vf->view_interval) PetscFunctionReturn(PETSC_SUCCESS);
1946:   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
1947:   PetscCall(VecView(tao->gradient, vf->viewer));
1948:   PetscCall(PetscViewerPopFormat(vf->viewer));
1949:   PetscFunctionReturn(PETSC_SUCCESS);
1950: }

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

1955:   Collective

1957:   Input Parameters:
1958: + tao - the `Tao` context
1959: - vf  - `PetscViewerAndFormat` context

1961:   Options Database Key:
1962: . -tao_monitor_step - view the step vector at each iteration

1964:   Level: advanced

1966: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1967: @*/
1968: PetscErrorCode TaoMonitorStep(Tao tao, PetscViewerAndFormat *vf)
1969: {
1970:   PetscFunctionBegin;
1972:   if (vf->view_interval > 0 && tao->niter % vf->view_interval) PetscFunctionReturn(PETSC_SUCCESS);
1973:   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
1974:   PetscCall(VecView(tao->stepdirection, vf->viewer));
1975:   PetscCall(PetscViewerPopFormat(vf->viewer));
1976:   PetscFunctionReturn(PETSC_SUCCESS);
1977: }

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

1982:   Collective

1984:   Input Parameters:
1985: + tao - the `Tao` context
1986: - ctx - `TaoMonitorDraw` context

1988:   Options Database Key:
1989: . -tao_monitor_solution_draw - draw the solution at each iteration

1991:   Level: advanced

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

1997: .seealso: [](ch_tao), `Tao`, `TaoMonitorSolution()`, `TaoMonitorSet()`, `TaoMonitorGradientDraw()`, `TaoMonitorDrawCtxCreate()`,
1998:           `TaoMonitorDrawCtxDestroy()`
1999: @*/
2000: PetscErrorCode TaoMonitorSolutionDraw(Tao tao, PetscCtx ctx)
2001: {
2002:   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;

2004:   PetscFunctionBegin;
2006:   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(PETSC_SUCCESS);
2007:   PetscCall(VecView(tao->solution, ictx->viewer));
2008:   PetscFunctionReturn(PETSC_SUCCESS);
2009: }

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

2014:   Collective

2016:   Input Parameters:
2017: + tao - the `Tao` context
2018: - ctx - `PetscViewer` context

2020:   Options Database Key:
2021: . -tao_monitor_gradient_draw - draw the gradient at each iteration

2023:   Level: advanced

2025: .seealso: [](ch_tao), `Tao`, `TaoMonitorGradient()`, `TaoMonitorSet()`, `TaoMonitorSolutionDraw()`
2026: @*/
2027: PetscErrorCode TaoMonitorGradientDraw(Tao tao, PetscCtx ctx)
2028: {
2029:   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;

2031:   PetscFunctionBegin;
2033:   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(PETSC_SUCCESS);
2034:   PetscCall(VecView(tao->gradient, ictx->viewer));
2035:   PetscFunctionReturn(PETSC_SUCCESS);
2036: }

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

2041:   Collective

2043:   Input Parameters:
2044: + tao - the `Tao` context
2045: - ctx - the `PetscViewer` context

2047:   Options Database Key:
2048: . -tao_monitor_step_draw - draw the step direction at each iteration

2050:   Level: advanced

2052: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorSolutionDraw`
2053: @*/
2054: PetscErrorCode TaoMonitorStepDraw(Tao tao, PetscCtx ctx)
2055: {
2056:   PetscViewer viewer = (PetscViewer)ctx;

2058:   PetscFunctionBegin;
2061:   PetscCall(VecView(tao->stepdirection, viewer));
2062:   PetscFunctionReturn(PETSC_SUCCESS);
2063: }

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

2068:   Collective

2070:   Input Parameters:
2071: + tao - the `Tao` context
2072: - vf  - `PetscViewerAndFormat` context

2074:   Options Database Key:
2075: . -tao_monitor_ls_residual - view the residual at each iteration

2077:   Level: advanced

2079: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
2080: @*/
2081: PetscErrorCode TaoMonitorResidual(Tao tao, PetscViewerAndFormat *vf)
2082: {
2083:   PetscFunctionBegin;
2085:   if (vf->view_interval > 0 && tao->niter % vf->view_interval) PetscFunctionReturn(PETSC_SUCCESS);
2086:   PetscCall(PetscViewerPushFormat(vf->viewer, vf->format));
2087:   PetscCall(VecView(tao->ls_res, vf->viewer));
2088:   PetscCall(PetscViewerPopFormat(vf->viewer));
2089:   PetscFunctionReturn(PETSC_SUCCESS);
2090: }

2092: /*@
2093:   TaoDefaultConvergenceTest - Determines whether the solver should continue iterating
2094:   or terminate.

2096:   Collective

2098:   Input Parameters:
2099: + tao   - the `Tao` context
2100: - dummy - unused dummy context

2102:   Level: developer

2104:   Notes:
2105:   This routine checks the residual in the optimality conditions, the
2106:   relative residual in the optimity conditions, the number of function
2107:   evaluations, and the function value to test convergence.  Some
2108:   solvers may use different convergence routines.

2110: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoGetConvergedReason()`, `TaoSetConvergedReason()`
2111: @*/
2112: PetscErrorCode TaoDefaultConvergenceTest(Tao tao, void *dummy)
2113: {
2114:   PetscInt           niter     = tao->niter, nfuncs;
2115:   PetscInt           max_funcs = tao->max_funcs;
2116:   PetscReal          gnorm = tao->residual, gnorm0 = tao->gnorm0;
2117:   PetscReal          f = tao->fc, steptol = tao->steptol, trradius = tao->step;
2118:   PetscReal          gatol = tao->gatol, grtol = tao->grtol, gttol = tao->gttol;
2119:   PetscReal          catol = tao->catol, crtol = tao->crtol;
2120:   PetscReal          fmin = tao->fmin, cnorm = tao->cnorm;
2121:   TaoConvergedReason reason = tao->reason;

2123:   PetscFunctionBegin;
2125:   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);

2127:   PetscCall(TaoGetCurrentFunctionEvaluations(tao, &nfuncs));
2128:   if (PetscIsInfOrNanReal(f)) {
2129:     PetscCall(PetscInfo(tao, "Failed to converged, function value is infinity or NaN\n"));
2130:     reason = TAO_DIVERGED_NAN;
2131:   } else if (f <= fmin && cnorm <= catol) {
2132:     PetscCall(PetscInfo(tao, "Converged due to function value %g < minimum function value %g\n", (double)f, (double)fmin));
2133:     reason = TAO_CONVERGED_MINF;
2134:   } else if (gnorm <= gatol && cnorm <= catol) {
2135:     PetscCall(PetscInfo(tao, "Converged due to residual norm ||g(X)||=%g < %g\n", (double)gnorm, (double)gatol));
2136:     reason = TAO_CONVERGED_GATOL;
2137:   } else if (f != 0 && PetscAbsReal(gnorm / f) <= grtol && cnorm <= crtol) {
2138:     PetscCall(PetscInfo(tao, "Converged due to residual ||g(X)||/|f(X)| =%g < %g\n", (double)(gnorm / f), (double)grtol));
2139:     reason = TAO_CONVERGED_GRTOL;
2140:   } else if (gnorm0 != 0 && ((gttol == 0 && gnorm == 0) || gnorm / gnorm0 < gttol) && cnorm <= crtol) {
2141:     PetscCall(PetscInfo(tao, "Converged due to relative residual norm ||g(X)||/||g(X0)|| = %g < %g\n", (double)(gnorm / gnorm0), (double)gttol));
2142:     reason = TAO_CONVERGED_GTTOL;
2143:   } else if (max_funcs != PETSC_UNLIMITED && nfuncs > max_funcs) {
2144:     PetscCall(PetscInfo(tao, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", nfuncs, max_funcs));
2145:     reason = TAO_DIVERGED_MAXFCN;
2146:   } else if (tao->lsflag != 0) {
2147:     PetscCall(PetscInfo(tao, "Tao Line Search failure.\n"));
2148:     reason = TAO_DIVERGED_LS_FAILURE;
2149:   } else if (trradius < steptol && niter > 0) {
2150:     PetscCall(PetscInfo(tao, "Trust region/step size too small: %g < %g\n", (double)trradius, (double)steptol));
2151:     reason = TAO_CONVERGED_STEPTOL;
2152:   } else if (niter >= tao->max_it) {
2153:     PetscCall(PetscInfo(tao, "Exceeded maximum number of iterations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", niter, tao->max_it));
2154:     reason = TAO_DIVERGED_MAXITS;
2155:   } else {
2156:     reason = TAO_CONTINUE_ITERATING;
2157:   }
2158:   tao->reason = reason;
2159:   PetscFunctionReturn(PETSC_SUCCESS);
2160: }

2162: /*@
2163:   TaoSetOptionsPrefix - Sets the prefix used for searching for all
2164:   Tao options in the database.

2166:   Logically Collective

2168:   Input Parameters:
2169: + tao - the `Tao` context
2170: - p   - the prefix string to prepend to all Tao option requests

2172:   Level: advanced

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

2178:   For example, to distinguish between the runtime options for two
2179:   different Tao solvers, one could call
2180: .vb
2181:       TaoSetOptionsPrefix(tao1,"sys1_")
2182:       TaoSetOptionsPrefix(tao2,"sys2_")
2183: .ve

2185:   This would enable use of different options for each system, such as
2186: .vb
2187:       -sys1_tao_method blmvm -sys1_tao_grtol 1.e-3
2188:       -sys2_tao_method lmvm  -sys2_tao_grtol 1.e-4
2189: .ve

2191: .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoAppendOptionsPrefix()`, `TaoGetOptionsPrefix()`
2192: @*/
2193: PetscErrorCode TaoSetOptionsPrefix(Tao tao, const char p[])
2194: {
2195:   PetscFunctionBegin;
2197:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tao, p));
2198:   if (tao->linesearch) PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, p));
2199:   if (tao->ksp) PetscCall(KSPSetOptionsPrefix(tao->ksp, p));
2200:   if (tao->callbacks) {
2201:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tao->callbacks, p));
2202:     PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao->callbacks, "callbacks_"));
2203:   }
2204:   PetscFunctionReturn(PETSC_SUCCESS);
2205: }

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

2210:   Logically Collective

2212:   Input Parameters:
2213: + tao - the `Tao` solver context
2214: - p   - the prefix string to prepend to all `Tao` option requests

2216:   Level: advanced

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

2222: .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoSetOptionsPrefix()`, `TaoGetOptionsPrefix()`
2223: @*/
2224: PetscErrorCode TaoAppendOptionsPrefix(Tao tao, const char p[])
2225: {
2226:   PetscFunctionBegin;
2228:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao, p));
2229:   if (tao->linesearch) PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao->linesearch, p));
2230:   if (tao->ksp) PetscCall(KSPAppendOptionsPrefix(tao->ksp, p));
2231:   if (tao->callbacks) {
2232:     const char *prefix;

2234:     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)tao, &prefix));
2235:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tao->callbacks, prefix));
2236:     PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao->callbacks, "callbacks_"));
2237:   }
2238:   PetscFunctionReturn(PETSC_SUCCESS);
2239: }

2241: /*@
2242:   TaoGetOptionsPrefix - Gets the prefix used for searching for all
2243:   Tao options in the database

2245:   Not Collective

2247:   Input Parameter:
2248: . tao - the `Tao` context

2250:   Output Parameter:
2251: . p - pointer to the prefix string used is returned

2253:   Level: advanced

2255: .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoSetOptionsPrefix()`, `TaoAppendOptionsPrefix()`
2256: @*/
2257: PetscErrorCode TaoGetOptionsPrefix(Tao tao, const char *p[])
2258: {
2259:   PetscFunctionBegin;
2261:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)tao, p));
2262:   PetscFunctionReturn(PETSC_SUCCESS);
2263: }

2265: /*@
2266:   TaoSetType - Sets the `TaoType` for the minimization solver.

2268:   Collective

2270:   Input Parameters:
2271: + tao  - the `Tao` solver context
2272: - type - a known method

2274:   Options Database Key:
2275: . -tao_type type - Sets the method; see `TaoType`

2277:   Level: intermediate

2279:   Note:
2280:   Calling this function resets the convergence test to `TaoDefaultConvergenceTest()`.
2281:   If a custom convergence test has been set with `TaoSetConvergenceTest()`, it must
2282:   be set again after calling `TaoSetType()`.

2284: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoGetType()`, `TaoType`
2285: @*/
2286: PetscErrorCode TaoSetType(Tao tao, TaoType type)
2287: {
2288:   PetscErrorCode (*create_xxx)(Tao);
2289:   PetscBool issame;

2291:   PetscFunctionBegin;

2294:   PetscCall(PetscObjectTypeCompare((PetscObject)tao, type, &issame));
2295:   if (issame) PetscFunctionReturn(PETSC_SUCCESS);

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

2300:   /* Destroy the existing solver information */
2301:   PetscTryTypeMethod(tao, destroy);
2302:   PetscCall(KSPDestroy(&tao->ksp));
2303:   PetscCall(TaoLineSearchDestroy(&tao->linesearch));

2305:   /* Reinitialize type-specific function pointers in TaoOps structure */
2306:   tao->ops->setup           = NULL;
2307:   tao->ops->computedual     = NULL;
2308:   tao->ops->solve           = NULL;
2309:   tao->ops->view            = NULL;
2310:   tao->ops->setfromoptions  = NULL;
2311:   tao->ops->destroy         = NULL;
2312:   tao->ops->convergencetest = TaoDefaultConvergenceTest;

2314:   tao->setupcalled           = PETSC_FALSE;
2315:   tao->uses_gradient         = PETSC_FALSE;
2316:   tao->uses_hessian_matrices = PETSC_FALSE;

2318:   PetscCall((*create_xxx)(tao));
2319:   PetscCall(PetscObjectChangeTypeName((PetscObject)tao, type));
2320:   PetscFunctionReturn(PETSC_SUCCESS);
2321: }

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

2326:   Not Collective, No Fortran Support

2328:   Input Parameters:
2329: + sname - name of a new user-defined solver
2330: - func  - routine to create `TaoType` specific method context

2332:   Calling sequence of `func`:
2333: . tao - the `Tao` object to be created

2335:   Example Usage:
2336: .vb
2337:    TaoRegister("my_solver", MySolverCreate);
2338: .ve

2340:   Then, your solver can be chosen with the procedural interface via
2341: .vb
2342:   TaoSetType(tao, "my_solver")
2343: .ve
2344:   or at runtime via the option
2345: .vb
2346:   -tao_type my_solver
2347: .ve

2349:   Level: advanced

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

2354: .seealso: [](ch_tao), `Tao`, `TaoSetType()`, `TaoRegisterAll()`, `TaoRegisterDestroy()`
2355: @*/
2356: PetscErrorCode TaoRegister(const char sname[], PetscErrorCode (*func)(Tao tao))
2357: {
2358:   PetscFunctionBegin;
2359:   PetscCall(TaoInitializePackage());
2360:   PetscCall(PetscFunctionListAdd(&TaoList, sname, func));
2361:   PetscFunctionReturn(PETSC_SUCCESS);
2362: }

2364: /*@C
2365:   TaoRegisterDestroy - Frees the list of minimization solvers that were
2366:   registered by `TaoRegister()`.

2368:   Not Collective

2370:   Level: advanced

2372: .seealso: [](ch_tao), `Tao`, `TaoRegisterAll()`, `TaoRegister()`
2373: @*/
2374: PetscErrorCode TaoRegisterDestroy(void)
2375: {
2376:   PetscFunctionBegin;
2377:   PetscCall(PetscFunctionListDestroy(&TaoList));
2378:   TaoRegisterAllCalled = PETSC_FALSE;
2379:   PetscFunctionReturn(PETSC_SUCCESS);
2380: }

2382: /*@
2383:   TaoGetIterationNumber - Gets the number of `TaoSolve()` iterations completed
2384:   at this time.

2386:   Not Collective

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

2391:   Output Parameter:
2392: . iter - iteration number

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

2397:   Level: intermediate

2399: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`, `TaoGetResidualNorm()`, `TaoGetObjective()`
2400: @*/
2401: PetscErrorCode TaoGetIterationNumber(Tao tao, PetscInt *iter)
2402: {
2403:   PetscFunctionBegin;
2405:   PetscAssertPointer(iter, 2);
2406:   *iter = tao->niter;
2407:   PetscFunctionReturn(PETSC_SUCCESS);
2408: }

2410: /*@
2411:   TaoGetResidualNorm - Gets the current value of the norm of the residual (gradient)
2412:   at this time.

2414:   Not Collective

2416:   Input Parameter:
2417: . tao - the `Tao` context

2419:   Output Parameter:
2420: . value - the current value

2422:   Level: intermediate

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

2428: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`, `TaoGetIterationNumber()`, `TaoGetObjective()`
2429: @*/
2430: PetscErrorCode TaoGetResidualNorm(Tao tao, PetscReal *value)
2431: {
2432:   PetscFunctionBegin;
2434:   PetscAssertPointer(value, 2);
2435:   *value = tao->residual;
2436:   PetscFunctionReturn(PETSC_SUCCESS);
2437: }

2439: /*@
2440:   TaoSetIterationNumber - Sets the current iteration number.

2442:   Logically Collective

2444:   Input Parameters:
2445: + tao  - the `Tao` context
2446: - iter - iteration number

2448:   Level: developer

2450: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
2451: @*/
2452: PetscErrorCode TaoSetIterationNumber(Tao tao, PetscInt iter)
2453: {
2454:   PetscFunctionBegin;
2457:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2458:   tao->niter = iter;
2459:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2460:   PetscFunctionReturn(PETSC_SUCCESS);
2461: }

2463: /*@
2464:   TaoGetTotalIterationNumber - Gets the total number of `TaoSolve()` iterations
2465:   completed. This number keeps accumulating if multiple solves
2466:   are called with the `Tao` object.

2468:   Not Collective

2470:   Input Parameter:
2471: . tao - the `Tao` context

2473:   Output Parameter:
2474: . iter - number of iterations

2476:   Level: intermediate

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

2482: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
2483: @*/
2484: PetscErrorCode TaoGetTotalIterationNumber(Tao tao, PetscInt *iter)
2485: {
2486:   PetscFunctionBegin;
2488:   PetscAssertPointer(iter, 2);
2489:   *iter = tao->ntotalits;
2490:   PetscFunctionReturn(PETSC_SUCCESS);
2491: }

2493: /*@
2494:   TaoSetTotalIterationNumber - Sets the current total iteration number.

2496:   Logically Collective

2498:   Input Parameters:
2499: + tao  - the `Tao` context
2500: - iter - the iteration number

2502:   Level: developer

2504: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
2505: @*/
2506: PetscErrorCode TaoSetTotalIterationNumber(Tao tao, PetscInt iter)
2507: {
2508:   PetscFunctionBegin;
2511:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2512:   tao->ntotalits = iter;
2513:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2514:   PetscFunctionReturn(PETSC_SUCCESS);
2515: }

2517: /*@
2518:   TaoSetConvergedReason - Sets the termination flag on a `Tao` object

2520:   Logically Collective

2522:   Input Parameters:
2523: + tao    - the `Tao` context
2524: - reason - the `TaoConvergedReason`

2526:   Level: intermediate

2528: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`
2529: @*/
2530: PetscErrorCode TaoSetConvergedReason(Tao tao, TaoConvergedReason reason)
2531: {
2532:   PetscFunctionBegin;
2535:   tao->reason = reason;
2536:   PetscFunctionReturn(PETSC_SUCCESS);
2537: }

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

2542:   Not Collective

2544:   Input Parameter:
2545: . tao - the `Tao` solver context

2547:   Output Parameter:
2548: . reason - value of `TaoConvergedReason`

2550:   Level: intermediate

2552: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetConvergenceTest()`, `TaoSetTolerances()`
2553: @*/
2554: PetscErrorCode TaoGetConvergedReason(Tao tao, TaoConvergedReason *reason)
2555: {
2556:   PetscFunctionBegin;
2558:   PetscAssertPointer(reason, 2);
2559:   *reason = tao->reason;
2560:   PetscFunctionReturn(PETSC_SUCCESS);
2561: }

2563: /*@
2564:   TaoGetSolutionStatus - Get the current iterate, objective value,
2565:   residual, infeasibility, and termination from a `Tao` object

2567:   Not Collective

2569:   Input Parameter:
2570: . tao - the `Tao` context

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

2580:   Level: intermediate

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

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

2587: .seealso: [](ch_tao), `TaoMonitor()`, `TaoGetConvergedReason()`
2588: @*/
2589: PetscErrorCode TaoGetSolutionStatus(Tao tao, PetscInt *its, PetscReal *f, PetscReal *gnorm, PetscReal *cnorm, PetscReal *xdiff, TaoConvergedReason *reason)
2590: {
2591:   PetscFunctionBegin;
2593:   if (its) *its = tao->niter;
2594:   if (f) *f = tao->fc;
2595:   if (gnorm) *gnorm = tao->residual;
2596:   if (cnorm) *cnorm = tao->cnorm;
2597:   if (reason) *reason = tao->reason;
2598:   if (xdiff) *xdiff = tao->step;
2599:   PetscFunctionReturn(PETSC_SUCCESS);
2600: }

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

2605:   Not Collective

2607:   Input Parameter:
2608: . tao - the `Tao` solver context

2610:   Output Parameter:
2611: . type - the `TaoType`

2613:   Level: intermediate

2615: .seealso: [](ch_tao), `Tao`, `TaoType`, `TaoSetType()`
2616: @*/
2617: PetscErrorCode TaoGetType(Tao tao, TaoType *type)
2618: {
2619:   PetscFunctionBegin;
2621:   PetscAssertPointer(type, 2);
2622:   *type = ((PetscObject)tao)->type_name;
2623:   PetscFunctionReturn(PETSC_SUCCESS);
2624: }

2626: /*@C
2627:   TaoMonitor - Monitor the solver and the current solution.  This
2628:   routine will record the iteration number and residual statistics,
2629:   and call any monitors specified by the user.

2631:   Input Parameters:
2632: + tao        - the `Tao` context
2633: . its        - the current iterate number (>=0)
2634: . f          - the current objective function value
2635: . res        - the gradient norm, square root of the duality gap, or other measure indicating distance from optimality.  This measure will be recorded and
2636:           used for some termination tests.
2637: . cnorm      - the infeasibility of the current solution with regard to the constraints.
2638: - steplength - multiple of the step direction added to the previous iterate.

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

2643:   Level: developer

2645: .seealso: [](ch_tao), `Tao`, `TaoGetConvergedReason()`, `TaoMonitorDefault()`, `TaoMonitorSet()`
2646: @*/
2647: PetscErrorCode TaoMonitor(Tao tao, PetscInt its, PetscReal f, PetscReal res, PetscReal cnorm, PetscReal steplength)
2648: {
2649:   PetscInt i;

2651:   PetscFunctionBegin;
2653:   tao->fc       = f;
2654:   tao->residual = res;
2655:   tao->cnorm    = cnorm;
2656:   tao->step     = steplength;
2657:   if (!its) {
2658:     tao->cnorm0 = cnorm;
2659:     tao->gnorm0 = res;
2660:   }
2661:   PetscCall(VecLockReadPush(tao->solution));
2662:   for (i = 0; i < tao->numbermonitors; i++) PetscCall((*tao->monitor[i])(tao, tao->monitorcontext[i]));
2663:   PetscCall(VecLockReadPop(tao->solution));
2664:   PetscFunctionReturn(PETSC_SUCCESS);
2665: }

2667: /*@
2668:   TaoSetConvergenceHistory - Sets the array used to hold the convergence history.

2670:   Logically Collective

2672:   Input Parameters:
2673: + tao   - the `Tao` solver context
2674: . obj   - array to hold objective value history
2675: . resid - array to hold residual history
2676: . cnorm - array to hold constraint violation history
2677: . lits  - integer array holds the number of linear iterations for each Tao iteration
2678: . na    - size of `obj`, `resid`, and `cnorm`
2679: - reset - `PETSC_TRUE` indicates each new minimization resets the history counter to zero,
2680:            else it continues storing new values for new minimizations after the old ones

2682:   Level: intermediate

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

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

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

2697: .seealso: [](ch_tao), `TaoGetConvergenceHistory()`
2698: @*/
2699: PetscErrorCode TaoSetConvergenceHistory(Tao tao, PetscReal obj[], PetscReal resid[], PetscReal cnorm[], PetscInt lits[], PetscInt na, PetscBool reset)
2700: {
2701:   PetscFunctionBegin;
2703:   if (obj) PetscAssertPointer(obj, 2);
2704:   if (resid) PetscAssertPointer(resid, 3);
2705:   if (cnorm) PetscAssertPointer(cnorm, 4);
2706:   if (lits) PetscAssertPointer(lits, 5);

2708:   if (na == PETSC_DECIDE || na == PETSC_CURRENT) na = 1000;
2709:   if (!obj && !resid && !cnorm && !lits) {
2710:     PetscCall(PetscCalloc4(na, &obj, na, &resid, na, &cnorm, na, &lits));
2711:     tao->hist_malloc = PETSC_TRUE;
2712:   }

2714:   tao->hist_obj   = obj;
2715:   tao->hist_resid = resid;
2716:   tao->hist_cnorm = cnorm;
2717:   tao->hist_lits  = lits;
2718:   tao->hist_max   = na;
2719:   tao->hist_reset = reset;
2720:   tao->hist_len   = 0;
2721:   PetscFunctionReturn(PETSC_SUCCESS);
2722: }

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

2727:   Collective

2729:   Input Parameter:
2730: . tao - the `Tao` context

2732:   Output Parameters:
2733: + obj   - array used to hold objective value history
2734: . resid - array used to hold residual history
2735: . cnorm - array used to hold constraint violation history
2736: . lits  - integer array used to hold linear solver iteration count
2737: - nhist - size of `obj`, `resid`, `cnorm`, and `lits`

2739:   Level: advanced

2741:   Notes:
2742:   This routine must be preceded by calls to `TaoSetConvergenceHistory()`
2743:   and `TaoSolve()`, otherwise it returns useless information.

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

2749:   Fortran Notes:
2750:   The calling sequence is
2751: .vb
2752:    call TaoGetConvergenceHistory(Tao tao, PetscInt nhist, PetscErrorCode ierr)
2753: .ve
2754:   In other words this gets the current number of entries in the history. Access the history through the array you passed to `TaoSetConvergenceHistory()`

2756: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetConvergenceHistory()`
2757: @*/
2758: PetscErrorCode TaoGetConvergenceHistory(Tao tao, PetscReal **obj, PetscReal **resid, PetscReal **cnorm, PetscInt **lits, PetscInt *nhist)
2759: {
2760:   PetscFunctionBegin;
2762:   if (obj) *obj = tao->hist_obj;
2763:   if (cnorm) *cnorm = tao->hist_cnorm;
2764:   if (resid) *resid = tao->hist_resid;
2765:   if (lits) *lits = tao->hist_lits;
2766:   if (nhist) *nhist = tao->hist_len;
2767:   PetscFunctionReturn(PETSC_SUCCESS);
2768: }

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

2774:   Logically Collective

2776:   Input Parameters:
2777: + tao - the `Tao` context
2778: - ctx - the user context

2780:   Level: intermediate

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

2787: .seealso: [](ch_tao), `Tao`, `TaoGetApplicationContext()`
2788: @*/
2789: PetscErrorCode TaoSetApplicationContext(Tao tao, PetscCtx ctx)
2790: {
2791:   PetscFunctionBegin;
2793:   tao->ctx = ctx;
2794:   PetscFunctionReturn(PETSC_SUCCESS);
2795: }

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

2800:   Not Collective

2802:   Input Parameter:
2803: . tao - the `Tao` context

2805:   Output Parameter:
2806: . ctx - a pointer to the user context

2808:   Level: intermediate

2810:   Fortran Note:
2811:   This only works when the context is a Fortran derived type or a `PetscObject`. Define `ctx` with
2812: .vb
2813:   type(tUsertype), pointer :: ctx
2814: .ve

2816: .seealso: [](ch_tao), `Tao`, `TaoSetApplicationContext()`
2817: @*/
2818: PetscErrorCode TaoGetApplicationContext(Tao tao, PetscCtxRt ctx)
2819: {
2820:   PetscFunctionBegin;
2822:   PetscAssertPointer(ctx, 2);
2823:   *(void **)ctx = tao->ctx;
2824:   PetscFunctionReturn(PETSC_SUCCESS);
2825: }

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

2830:   Collective

2832:   Input Parameters:
2833: + tao - the `Tao` context
2834: - M   - matrix that defines the norm

2836:   Level: beginner

2838: .seealso: [](ch_tao), `Tao`, `TaoGetGradientNorm()`, `TaoGradientNorm()`
2839: @*/
2840: PetscErrorCode TaoSetGradientNorm(Tao tao, Mat M)
2841: {
2842:   PetscFunctionBegin;
2845:   PetscCall(PetscObjectReference((PetscObject)M));
2846:   PetscCall(MatDestroy(&tao->gradient_norm));
2847:   PetscCall(VecDestroy(&tao->gradient_norm_tmp));
2848:   tao->gradient_norm = M;
2849:   PetscCall(MatCreateVecs(M, NULL, &tao->gradient_norm_tmp));
2850:   PetscFunctionReturn(PETSC_SUCCESS);
2851: }

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

2856:   Not Collective

2858:   Input Parameter:
2859: . tao - the `Tao` context

2861:   Output Parameter:
2862: . M - gradient norm

2864:   Level: beginner

2866: .seealso: [](ch_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGradientNorm()`
2867: @*/
2868: PetscErrorCode TaoGetGradientNorm(Tao tao, Mat *M)
2869: {
2870:   PetscFunctionBegin;
2872:   PetscAssertPointer(M, 2);
2873:   *M = tao->gradient_norm;
2874:   PetscFunctionReturn(PETSC_SUCCESS);
2875: }

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

2880:   Collective

2882:   Input Parameters:
2883: + tao      - the `Tao` context
2884: . gradient - the gradient
2885: - type     - the norm type

2887:   Output Parameter:
2888: . gnorm - the gradient norm

2890:   Level: advanced

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

2895:   Developer Notes:
2896:   Should be named `TaoComputeGradientNorm()`.

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

2901: .seealso: [](ch_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGetGradientNorm()`
2902: @*/
2903: PetscErrorCode TaoGradientNorm(Tao tao, Vec gradient, NormType type, PetscReal *gnorm)
2904: {
2905:   PetscFunctionBegin;
2909:   PetscAssertPointer(gnorm, 4);
2910:   if (tao->gradient_norm) {
2911:     PetscScalar gnorms;

2913:     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.");
2914:     PetscCall(MatMult(tao->gradient_norm, gradient, tao->gradient_norm_tmp));
2915:     PetscCall(VecDot(gradient, tao->gradient_norm_tmp, &gnorms));
2916:     *gnorm = PetscRealPart(PetscSqrtScalar(gnorms));
2917:   } else {
2918:     PetscCall(VecNorm(gradient, type, gnorm));
2919:   }
2920:   PetscFunctionReturn(PETSC_SUCCESS);
2921: }

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

2926:   Collective

2928:   Input Parameters:
2929: + comm     - the communicator to share the context
2930: . host     - the name of the X Windows host that will display the monitor
2931: . label    - the label to put at the top of the display window
2932: . x        - the horizontal coordinate of the lower left corner of the window to open
2933: . y        - the vertical coordinate of the lower left corner of the window to open
2934: . m        - the width of the window
2935: . n        - the height of the window
2936: - howoften - how many `Tao` iterations between displaying the monitor information

2938:   Output Parameter:
2939: . ctx - the monitor context

2941:   Options Database Keys:
2942: + -tao_monitor_solution_draw - use `TaoMonitorSolutionDraw()` to monitor the solution
2943: - -tao_draw_solution_initial - show initial guess as well as current solution

2945:   Level: intermediate

2947:   Note:
2948:   The context this creates, along with `TaoMonitorSolutionDraw()`, and `TaoMonitorDrawCtxDestroy()`
2949:   are passed to `TaoMonitorSet()`.

2951: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorDrawCtx()`
2952: @*/
2953: PetscErrorCode TaoMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TaoMonitorDrawCtx *ctx)
2954: {
2955:   PetscFunctionBegin;
2956:   PetscCall(PetscNew(ctx));
2957:   PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
2958:   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
2959:   (*ctx)->howoften = howoften;
2960:   PetscFunctionReturn(PETSC_SUCCESS);
2961: }

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

2966:   Collective

2968:   Input Parameter:
2969: . ictx - the monitor context

2971:   Level: intermediate

2973:   Note:
2974:   This is passed to `TaoMonitorSet()` as the final argument, along with `TaoMonitorSolutionDraw()`, and the context
2975:   obtained with `TaoMonitorDrawCtxCreate()`.

2977: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorSolutionDraw()`
2978: @*/
2979: PetscErrorCode TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx *ictx)
2980: {
2981:   PetscFunctionBegin;
2982:   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
2983:   PetscCall(PetscFree(*ictx));
2984:   PetscFunctionReturn(PETSC_SUCCESS);
2985: }

2987: /*@
2988:   TaoGetTerm - Get the entire objective function of the `Tao` as a
2989:   single `TaoTerm` in the form $\alpha f(Ax; p)$, where $\alpha$ is a scaling
2990:   coefficient, $f$ is a `TaoTerm`, $A$ is an (optional) map and $p$ are the parameters of $f$.

2992:   Not collective

2994:   Input Parameter:
2995: . tao - a `Tao` context

2997:   Output Parameters:
2998: + scale  - the scale of the term
2999: . term   - a `TaoTerm` for the real-valued function defining the objective
3000: . params - the vector of parameters for `term`, or `NULL` if no parameters were specified for `term`
3001: - map    - a map from the solution space of `tao` to the solution space of `term`, if `NULL` then the map is the identity

3003:   Level: intermediate

3005:   Notes:
3006:   If the objective function was defined by providing function callbacks directly to `Tao` (for example, with `TaoSetObjectiveAndGradient()`), then
3007:   `TaoGetTerm` will return a `TaoTerm` with the type `TAOTERMCALLBACKS` that encapsulates
3008:   those functions.

3010:   If multiple `TaoTerms` were provided to `Tao` via, for example, `TaoAddTerm()`, or in combination with giving functions directly to `Tao`, then the type `TAOTERMSUM` is returned.

3012: .seealso: [](ch_tao), `Tao`, `TaoTerm`, `TAOTERMSUM`, `TaoAddTerm()`
3013: @*/
3014: PetscErrorCode TaoGetTerm(Tao tao, PetscReal *scale, TaoTerm *term, Vec *params, Mat *map)
3015: {
3016:   PetscFunctionBegin;
3018:   if (scale) PetscAssertPointer(scale, 2);
3019:   if (term) PetscAssertPointer(term, 3);
3020:   if (params) PetscAssertPointer(params, 4);
3021:   if (map) PetscAssertPointer(map, 5);
3022:   PetscCall(TaoTermMappingGetData(&tao->objective_term, NULL, scale, term, map));
3023:   if (params) *params = tao->objective_parameters;
3024:   PetscFunctionReturn(PETSC_SUCCESS);
3025: }

3027: /*@
3028:   TaoAddTerm - Add a `term` to the objective function. If `Tao` is empty,
3029:   `term` will be the objective of `Tao`.

3031:   Collective

3033:   Input Parameters:
3034: + tao    - a `Tao` solver context
3035: . prefix - the prefix used for configuring the new term (if `NULL`, the index of the term will be used as a prefix, e.g. "0_", "1_", etc.)
3036: . scale  - scaling coefficient for the new term
3037: . term   - the real-valued function defining the new term
3038: . params - (optional) parameters for the new term.  It is up to each implementation of `TaoTerm` to determine how it behaves when parameters are omitted.
3039: - map    - (optional) a map from the `tao` solution space to the `term` solution space; if `NULL` the map is assumed to be the identity

3041:   Level: beginner

3043:   Notes:
3044:   If the objective function was $f(x)$, after calling `TaoAddTerm()` it becomes
3045:   $f(x) + \alpha g(Ax; p)$, where $\alpha$ is the `scale`, $g$ is the `term`, $A$ is the
3046:   (optional) `map`, and $p$ are the (optional) `params` of $g$.

3048:   The `map` $A$ transforms the `Tao` solution vector into the term's solution space.
3049:   For example, if the `Tao` solution vector is $x \in \mathbb{R}^n$ and the mapping
3050:   matrix is $A \in \mathbb{R}^{m \times n}$, then the term evaluates $g(Ax; p)$ with
3051:   $Ax \in \mathbb{R}^m$. The term's solution space is therefore $\mathbb{R}^m$. If the map is
3052:   `NULL`, the identity is used and the term's solution space must match the `Tao` solution space.
3053:   `Tao` automatically applies the chain rule for gradients ($A^T \nabla g$) and Hessians
3054:   ($A^T \nabla^2 g \, A$) with respect to $x$.

3056:   The `params` $p$ are fixed data that are not optimized over. Some `TaoTermType`s
3057:   require the parameter space to be related to the term's solution space (e.g., the same
3058:   size); when a mapping matrix $A$ is used, the parameter space may depend on either the row
3059:   or column space of $A$.  See the documentation for each `TaoTermType`.

3061:   Currently, `TaoAddTerm()` does not support bounded Newton solvers (`TAOBNK`,`TAOBNLS`,`TAOBNTL`,`TAOBNTR`,and `TAOBQNK`)

3063: .seealso: [](ch_tao), `Tao`, `TaoTerm`, `TAOTERMSUM`, `TaoGetTerm()`
3064: @*/
3065: PetscErrorCode TaoAddTerm(Tao tao, const char prefix[], PetscReal scale, TaoTerm term, Vec params, Mat map)
3066: {
3067:   PetscBool is_sum, is_callback;
3068:   PetscInt  num_old_terms;
3069:   Vec      *vec_list = NULL;

3071:   PetscFunctionBegin;
3073:   if (prefix) PetscAssertPointer(prefix, 2);
3076:   PetscCheckSameComm(tao, 1, term, 4);
3077:   if (params) {
3079:     PetscCheckSameComm(tao, 1, params, 5);
3080:   }
3081:   if (map) {
3083:     PetscCheckSameComm(tao, 1, map, 6);
3084:   }
3085:   // If user is using TaoAddTerm, before setting any terms or callbacks,
3086:   // then tao->objective_term.term is empty callback, which we want to remove.
3087:   PetscCall(PetscObjectTypeCompare((PetscObject)tao->objective_term.term, TAOTERMCALLBACKS, &is_callback));
3088:   PetscCall(PetscObjectTypeCompare((PetscObject)term, TAOTERMSUM, &is_sum));
3089:   PetscCheck(!is_sum, PetscObjectComm((PetscObject)term), PETSC_ERR_ARG_WRONG, "TaoAddTerm does not support adding TAOTERMSUM");
3090:   if (is_callback) {
3091:     PetscBool is_obj, is_objgrad, is_grad;

3093:     PetscCall(TaoTermIsObjectiveDefined(tao->objective_term.term, &is_obj));
3094:     PetscCall(TaoTermIsObjectiveAndGradientDefined(tao->objective_term.term, &is_objgrad));
3095:     PetscCall(TaoTermIsGradientDefined(tao->objective_term.term, &is_grad));
3096:     // Empty callback term
3097:     if (!(is_obj || is_objgrad || is_grad)) {
3098:       PetscCall(TaoTermMappingSetData(&tao->objective_term, NULL, scale, term, map));
3099:       PetscCall(PetscObjectReference((PetscObject)params));
3100:       PetscCall(VecDestroy(&tao->objective_parameters));
3101:       // Empty callback term. Destroy hessians, as they are not needed
3102:       PetscCall(MatDestroy(&tao->hessian));
3103:       PetscCall(MatDestroy(&tao->hessian_pre));
3104:       tao->objective_parameters = params;
3105:       tao->term_set             = PETSC_TRUE;
3106:       PetscFunctionReturn(PETSC_SUCCESS);
3107:     }
3108:   }
3109:   PetscCall(PetscObjectTypeCompare((PetscObject)tao->objective_term.term, TAOTERMSUM, &is_sum));
3110:   // One TaoTerm has been set. Create TAOTERMSUM to store that, and the new one
3111:   if (!is_sum) {
3112:     TaoTerm     old_sum;
3113:     const char *tao_prefix;
3114:     const char *term_prefix;

3116:     PetscCall(TaoTermDuplicate(tao->objective_term.term, TAOTERM_DUPLICATE_SIZEONLY, &old_sum));
3117:     if (tao->objective_term.map) {
3118:       VecType     map_vectype;
3119:       VecType     param_vectype;
3120:       PetscLayout cmap, param_layout;

3122:       PetscCall(MatGetVecType(tao->objective_term.map, &map_vectype));
3123:       PetscCall(MatGetLayouts(tao->objective_term.map, NULL, &cmap));
3124:       PetscCall(TaoTermGetParametersVecType(old_sum, &param_vectype));
3125:       PetscCall(TaoTermGetParametersLayout(old_sum, &param_layout));

3127:       PetscCall(TaoTermSetSolutionVecType(old_sum, map_vectype));
3128:       PetscCall(TaoTermSetParametersVecType(old_sum, param_vectype));
3129:       PetscCall(TaoTermSetSolutionLayout(old_sum, cmap));
3130:       PetscCall(TaoTermSetParametersLayout(old_sum, param_layout));
3131:     }

3133:     PetscCall(TaoTermSetType(old_sum, TAOTERMSUM));
3134:     PetscCall(TaoGetOptionsPrefix(tao, &tao_prefix));
3135:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)old_sum, tao_prefix));
3136:     PetscCall(TaoTermSumSetNumberTerms(old_sum, 1));
3137:     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)tao->objective_term.term, &term_prefix));
3138:     PetscCall(TaoTermSumSetTerm(old_sum, 0, term_prefix, tao->objective_term.scale, tao->objective_term.term, tao->objective_term.map));
3139:     PetscCall(TaoTermSumSetTermHessianMatrices(old_sum, 0, NULL, NULL, tao->hessian, tao->hessian_pre));
3140:     PetscCall(MatDestroy(&tao->hessian));
3141:     PetscCall(MatDestroy(&tao->hessian_pre));
3142:     PetscCall(TaoTermMappingReset(&tao->objective_term));
3143:     PetscCall(TaoTermMappingSetData(&tao->objective_term, NULL, 1.0, old_sum, NULL));
3144:     if (tao->objective_parameters) {
3145:       // convert the parameters to a VECNEST
3146:       Vec subvecs[1];

3148:       subvecs[0]                = tao->objective_parameters;
3149:       tao->objective_parameters = NULL;
3150:       PetscCall(TaoTermSumParametersPack(old_sum, subvecs, &tao->objective_parameters));
3151:       PetscCall(VecDestroy(&subvecs[0]));
3152:     }
3153:     PetscCall(TaoTermDestroy(&old_sum));
3154:     tao->num_terms = 1;
3155:   }
3156:   PetscCall(TaoTermSumGetNumberTerms(tao->objective_term.term, &num_old_terms));
3157:   if (tao->objective_parameters || params) {
3158:     PetscCall(PetscCalloc1(num_old_terms + 1, &vec_list));
3159:     if (tao->objective_parameters) PetscCall(TaoTermSumParametersUnpack(tao->objective_term.term, &tao->objective_parameters, vec_list));
3160:     PetscCall(PetscObjectReference((PetscObject)params));
3161:     vec_list[num_old_terms] = params;
3162:   }
3163:   PetscCall(TaoTermSumAddTerm(tao->objective_term.term, prefix, scale, term, map, NULL));
3164:   tao->num_terms++;
3165:   if (vec_list) {
3166:     PetscInt num_terms = num_old_terms + 1;
3167:     PetscCall(TaoTermSumParametersPack(tao->objective_term.term, vec_list, &tao->objective_parameters));
3168:     for (PetscInt i = 0; i < num_terms; i++) PetscCall(VecDestroy(&vec_list[i]));
3169:     PetscCall(PetscFree(vec_list));
3170:   }
3171:   PetscFunctionReturn(PETSC_SUCCESS);
3172: }