Actual source code: taosolver.c

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

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

  8: PetscClassId TAO_CLASSID = 0;

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

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

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

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

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

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

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

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

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

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

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

 74:   Collective

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

 79:   Level: developer

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

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

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

105:   Collective

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

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

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

116:   Level: beginner

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

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

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

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

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

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

146:   Collective

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

151:   Level: beginner

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

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

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

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

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

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

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

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

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

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

215:   Collective

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

220:   Level: advanced

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

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

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

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

261:   Collective

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

266:   Level: beginner

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

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

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

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

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

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

343:   Logically Collective

345:   Input Parameters:
346: + tao  - Tao context
347: - flag - `PETSC_TRUE` or `PETSC_FALSE`

349:   Level: advanced

351:   Note:
352:   See `SNESKSPSetUseEW()` for customization details.

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

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

368:   Collective

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

377:   Level: developer

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

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

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

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

410: /*@
411:   TaoSetFromOptions - Sets various Tao parameters from the options database

413:   Collective

415:   Input Parameter:
416: . tao - the `Tao` solver context

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

448:   Level: beginner

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

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

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

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

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

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

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

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

492:   gatol = tao->gatol;
493:   grtol = tao->grtol;
494:   gttol = tao->gttol;
495:   PetscCall(PetscOptionsReal("-tao_gatol", "Stop if norm of gradient less than", "TaoSetTolerances", tao->gatol, &gatol, NULL));
496:   PetscCall(PetscOptionsReal("-tao_grtol", "Stop if norm of gradient divided by the function value is less than", "TaoSetTolerances", tao->grtol, &grtol, NULL));
497:   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));
498:   PetscCall(TaoSetTolerances(tao, gatol, grtol, gttol));

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

627:   PetscTryTypeMethod(tao, setfromoptions, PetscOptionsObject);

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

633:   if (tao->linesearch) PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
634:   PetscFunctionReturn(PETSC_SUCCESS);
635: }

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

640:   Collective

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

647:   Options Database Key:
648: . -name [viewertype][:...] - option name and values. See `PetscObjectViewFromOptions()` for the possible arguments

650:   Level: intermediate

652: .seealso: [](ch_tao), `Tao`, `TaoView`, `PetscObjectViewFromOptions()`, `TaoCreate()`
653: @*/
654: PetscErrorCode TaoViewFromOptions(Tao A, PetscObject obj, const char name[])
655: {
656:   PetscFunctionBegin;
658:   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
659:   PetscFunctionReturn(PETSC_SUCCESS);
660: }

662: /*@
663:   TaoView - Prints information about the `Tao` object

665:   Collective

667:   Input Parameters:
668: + tao    - the `Tao` context
669: - viewer - visualization context

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

674:   Level: beginner

676:   Notes:
677:   The available visualization contexts include
678: +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
679: -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
680:   output where only the first processor opens
681:   the file.  All other processors send their
682:   data to the first processor to print.

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

687: .seealso: [](ch_tao), `Tao`, `PetscViewerASCIIOpen()`
688: @*/
689: PetscErrorCode TaoView(Tao tao, PetscViewer viewer)
690: {
691:   PetscBool isascii, isstring;
692:   TaoType   type;

694:   PetscFunctionBegin;
696:   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(((PetscObject)tao)->comm, &viewer));
698:   PetscCheckSameComm(tao, 1, viewer, 2);

700:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
701:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
702:   if (isascii) {
703:     PetscViewerFormat format;

705:     PetscCall(PetscViewerGetFormat(viewer, &format));
706:     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)tao, viewer));

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

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

758:     PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: gatol=%g,", (double)tao->gatol));
759:     PetscCall(PetscViewerASCIIPrintf(viewer, " grtol=%g,", (double)tao->grtol));
760:     PetscCall(PetscViewerASCIIPrintf(viewer, " steptol=%g,", (double)tao->steptol));
761:     PetscCall(PetscViewerASCIIPrintf(viewer, " gttol=%g\n", (double)tao->gttol));
762:     PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Function/Gradient:=%g\n", (double)tao->residual));

764:     if (tao->constrained) {
765:       PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances:"));
766:       PetscCall(PetscViewerASCIIPrintf(viewer, " catol=%g,", (double)tao->catol));
767:       PetscCall(PetscViewerASCIIPrintf(viewer, " crtol=%g\n", (double)tao->crtol));
768:       PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Constraints:=%g\n", (double)tao->cnorm));
769:     }

771:     if (tao->trust < tao->steptol) {
772:       PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: steptol=%g\n", (double)tao->steptol));
773:       PetscCall(PetscViewerASCIIPrintf(viewer, "Final trust region radius:=%g\n", (double)tao->trust));
774:     }

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

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

782:     if (tao->objective_term.term->nobj > 0) {
783:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function evaluations=%" PetscInt_FMT ",", tao->objective_term.term->nobj));
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->ngrad > 0) {
788:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of gradient evaluations=%" PetscInt_FMT ",", tao->objective_term.term->ngrad));
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->objective_term.term->nobjgrad > 0) {
793:       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function/gradient evaluations=%" PetscInt_FMT ",", tao->objective_term.term->nobjgrad));
794:       if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, "    (max: unlimited)\n"));
795:       else PetscCall(PetscViewerASCIIPrintf(viewer, "    (max: %" PetscInt_FMT ")\n", tao->max_funcs));
796:     }
797:     if (tao->nres > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of residual evaluations=%" PetscInt_FMT "\n", tao->nres));
798:     if (tao->objective_term.term->nhess > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Hessian evaluations=%" PetscInt_FMT "\n", tao->objective_term.term->nhess));
799:     if (tao->nconstraints > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of constraint function evaluations=%" PetscInt_FMT "\n", tao->nconstraints));
800:     if (tao->njac > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Jacobian evaluations=%" PetscInt_FMT "\n", tao->njac));

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

863: /*@
864:   TaoSetRecycleHistory - Sets the boolean flag to enable/disable re-using
865:   iterate information from the previous `TaoSolve()`. This feature is disabled by
866:   default.

868:   Logically Collective

870:   Input Parameters:
871: + tao     - the `Tao` context
872: - recycle - boolean flag

874:   Options Database Key:
875: . -tao_recycle_history (true|false) - reuse the history

877:   Level: intermediate

879:   Notes:
880:   For conjugate gradient methods (`TAOBNCG`), this re-uses the latest search direction
881:   from the previous `TaoSolve()` call when computing the first search direction in a
882:   new solution. By default, CG methods set the first search direction to the
883:   negative gradient.

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

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

892: .seealso: [](ch_tao), `Tao`, `TaoGetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
893: @*/
894: PetscErrorCode TaoSetRecycleHistory(Tao tao, PetscBool recycle)
895: {
896:   PetscFunctionBegin;
899:   tao->recycle = recycle;
900:   PetscFunctionReturn(PETSC_SUCCESS);
901: }

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

907:   Logically Collective

909:   Input Parameter:
910: . tao - the `Tao` context

912:   Output Parameter:
913: . recycle - boolean flag

915:   Level: intermediate

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

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

931:   Logically Collective

933:   Input Parameters:
934: + tao   - the `Tao` context
935: . gatol - stop if norm of gradient is less than this
936: . grtol - stop if relative norm of gradient is less than this
937: - gttol - stop if norm of gradient is reduced by this factor

939:   Options Database Keys:
940: + -tao_gatol gatol - Sets gatol
941: . -tao_grtol grtol - Sets grtol
942: - -tao_gttol gttol - Sets gttol

944:   Stopping Criteria\:
945: .vb
946:   ||g(X)||                            <= gatol
947:   ||g(X)|| / |f(X)|                   <= grtol
948:   ||g(X)|| / ||g(X0)||                <= gttol
949: .ve

951:   Level: beginner

953:   Notes:
954:   Use `PETSC_CURRENT` to leave one or more tolerances unchanged.

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

958:   Fortran Note:
959:   Use `PETSC_CURRENT_REAL` or `PETSC_DETERMINE_REAL`

961: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`
962: @*/
963: PetscErrorCode TaoSetTolerances(Tao tao, PetscReal gatol, PetscReal grtol, PetscReal gttol)
964: {
965:   PetscFunctionBegin;

971:   if (gatol == (PetscReal)PETSC_DETERMINE) {
972:     tao->gatol = tao->default_gatol;
973:   } else if (gatol != (PetscReal)PETSC_CURRENT) {
974:     PetscCheck(gatol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative gatol not allowed");
975:     tao->gatol = gatol;
976:   }

978:   if (grtol == (PetscReal)PETSC_DETERMINE) {
979:     tao->grtol = tao->default_grtol;
980:   } else if (grtol != (PetscReal)PETSC_CURRENT) {
981:     PetscCheck(grtol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative grtol not allowed");
982:     tao->grtol = grtol;
983:   }

985:   if (gttol == (PetscReal)PETSC_DETERMINE) {
986:     tao->gttol = tao->default_gttol;
987:   } else if (gttol != (PetscReal)PETSC_CURRENT) {
988:     PetscCheck(gttol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative gttol not allowed");
989:     tao->gttol = gttol;
990:   }
991:   PetscFunctionReturn(PETSC_SUCCESS);
992: }

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

997:   Logically Collective

999:   Input Parameters:
1000: + tao   - the `Tao` context
1001: . catol - absolute constraint tolerance, constraint norm must be less than `catol` for used for `gatol` convergence criteria
1002: - crtol - relative constraint tolerance, constraint norm must be less than `crtol` for used for `gatol`, `gttol` convergence criteria

1004:   Options Database Keys:
1005: + -tao_catol catol - Sets catol
1006: - -tao_crtol crtol - Sets crtol

1008:   Level: intermediate

1010:   Notes:
1011:   Use `PETSC_CURRENT` to leave one or tolerance unchanged.

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

1015:   Fortran Note:
1016:   Use `PETSC_CURRENT_REAL` or `PETSC_DETERMINE_REAL`

1018: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`, `TaoGetConstraintTolerances()`, `TaoSetTolerances()`
1019: @*/
1020: PetscErrorCode TaoSetConstraintTolerances(Tao tao, PetscReal catol, PetscReal crtol)
1021: {
1022:   PetscFunctionBegin;

1027:   if (catol == (PetscReal)PETSC_DETERMINE) {
1028:     tao->catol = tao->default_catol;
1029:   } else if (catol != (PetscReal)PETSC_CURRENT) {
1030:     PetscCheck(catol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative catol not allowed");
1031:     tao->catol = catol;
1032:   }

1034:   if (crtol == (PetscReal)PETSC_DETERMINE) {
1035:     tao->crtol = tao->default_crtol;
1036:   } else if (crtol != (PetscReal)PETSC_CURRENT) {
1037:     PetscCheck(crtol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative crtol not allowed");
1038:     tao->crtol = crtol;
1039:   }
1040:   PetscFunctionReturn(PETSC_SUCCESS);
1041: }

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

1046:   Not Collective

1048:   Input Parameter:
1049: . tao - the `Tao` context

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

1055:   Level: intermediate

1057: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`, `TaoSetTolerances()`, `TaoSetConstraintTolerances()`
1058: @*/
1059: PetscErrorCode TaoGetConstraintTolerances(Tao tao, PetscReal *catol, PetscReal *crtol)
1060: {
1061:   PetscFunctionBegin;
1063:   if (catol) *catol = tao->catol;
1064:   if (crtol) *crtol = tao->crtol;
1065:   PetscFunctionReturn(PETSC_SUCCESS);
1066: }

1068: /*@
1069:   TaoSetFunctionLowerBound - Sets a bound on the solution objective value.
1070:   When an approximate solution with an objective value below this number
1071:   has been found, the solver will terminate.

1073:   Logically Collective

1075:   Input Parameters:
1076: + tao  - the Tao solver context
1077: - fmin - the tolerance

1079:   Options Database Key:
1080: . -tao_fmin fmin - sets the minimum function value

1082:   Level: intermediate

1084: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetTolerances()`
1085: @*/
1086: PetscErrorCode TaoSetFunctionLowerBound(Tao tao, PetscReal fmin)
1087: {
1088:   PetscFunctionBegin;
1091:   tao->fmin = fmin;
1092:   PetscFunctionReturn(PETSC_SUCCESS);
1093: }

1095: /*@
1096:   TaoGetFunctionLowerBound - Gets the bound on the solution objective value.
1097:   When an approximate solution with an objective value below this number
1098:   has been found, the solver will terminate.

1100:   Not Collective

1102:   Input Parameter:
1103: . tao - the `Tao` solver context

1105:   Output Parameter:
1106: . fmin - the minimum function value

1108:   Level: intermediate

1110: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetFunctionLowerBound()`
1111: @*/
1112: PetscErrorCode TaoGetFunctionLowerBound(Tao tao, PetscReal *fmin)
1113: {
1114:   PetscFunctionBegin;
1116:   PetscAssertPointer(fmin, 2);
1117:   *fmin = tao->fmin;
1118:   PetscFunctionReturn(PETSC_SUCCESS);
1119: }

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

1124:   Logically Collective

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

1130:   Options Database Key:
1131: . -tao_max_funcs nfcn - sets the maximum number of function evaluations

1133:   Level: intermediate

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

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

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

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

1162:   Logically Collective

1164:   Input Parameter:
1165: . tao - the `Tao` solver context

1167:   Output Parameter:
1168: . nfcn - the maximum number of function evaluations

1170:   Level: intermediate

1172: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1173: @*/
1174: PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao tao, PetscInt *nfcn)
1175: {
1176:   PetscFunctionBegin;
1178:   PetscAssertPointer(nfcn, 2);
1179:   *nfcn = tao->max_funcs;
1180:   PetscFunctionReturn(PETSC_SUCCESS);
1181: }

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

1186:   Not Collective

1188:   Input Parameter:
1189: . tao - the `Tao` solver context

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

1194:   Level: intermediate

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

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

1210:   Logically Collective

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

1216:   Options Database Key:
1217: . -tao_max_it its - sets the maximum number of iterations

1219:   Level: intermediate

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

1224:   Developer Note:
1225:   Also accepts the deprecated negative values to indicate no limit

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

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

1248:   Not Collective

1250:   Input Parameter:
1251: . tao - the `Tao` solver context

1253:   Output Parameter:
1254: . maxits - the maximum number of iterates

1256:   Level: intermediate

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

1269: /*@
1270:   TaoSetInitialTrustRegionRadius - Sets the initial trust region radius.

1272:   Logically Collective

1274:   Input Parameters:
1275: + tao    - a `Tao` optimization solver
1276: - radius - the trust region radius

1278:   Options Database Key:
1279: . -tao_trust0 radius - sets initial trust region radius

1281:   Level: intermediate

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

1286: .seealso: [](ch_tao), `Tao`, `TaoGetTrustRegionRadius()`, `TaoSetTrustRegionTolerance()`, `TAONTR`
1287: @*/
1288: PetscErrorCode TaoSetInitialTrustRegionRadius(Tao tao, PetscReal radius)
1289: {
1290:   PetscFunctionBegin;
1293:   if (radius == PETSC_DETERMINE) {
1294:     tao->trust0 = tao->default_trust0;
1295:   } else {
1296:     PetscCheck(radius > 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Radius must be positive");
1297:     tao->trust0 = radius;
1298:   }
1299:   PetscFunctionReturn(PETSC_SUCCESS);
1300: }

1302: /*@
1303:   TaoGetInitialTrustRegionRadius - Gets the initial trust region radius.

1305:   Not Collective

1307:   Input Parameter:
1308: . tao - a `Tao` optimization solver

1310:   Output Parameter:
1311: . radius - the trust region radius

1313:   Level: intermediate

1315: .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetCurrentTrustRegionRadius()`, `TAONTR`
1316: @*/
1317: PetscErrorCode TaoGetInitialTrustRegionRadius(Tao tao, PetscReal *radius)
1318: {
1319:   PetscFunctionBegin;
1321:   PetscAssertPointer(radius, 2);
1322:   *radius = tao->trust0;
1323:   PetscFunctionReturn(PETSC_SUCCESS);
1324: }

1326: /*@
1327:   TaoGetCurrentTrustRegionRadius - Gets the current trust region radius.

1329:   Not Collective

1331:   Input Parameter:
1332: . tao - a `Tao` optimization solver

1334:   Output Parameter:
1335: . radius - the trust region radius

1337:   Level: intermediate

1339: .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetInitialTrustRegionRadius()`, `TAONTR`
1340: @*/
1341: PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao tao, PetscReal *radius)
1342: {
1343:   PetscFunctionBegin;
1345:   PetscAssertPointer(radius, 2);
1346:   *radius = tao->trust;
1347:   PetscFunctionReturn(PETSC_SUCCESS);
1348: }

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

1353:   Not Collective

1355:   Input Parameter:
1356: . tao - the `Tao` context

1358:   Output Parameters:
1359: + gatol - stop if norm of gradient is less than this
1360: . grtol - stop if relative norm of gradient is less than this
1361: - gttol - stop if norm of gradient is reduced by a this factor

1363:   Level: intermediate

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

1368: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`
1369: @*/
1370: PetscErrorCode TaoGetTolerances(Tao tao, PetscReal *gatol, PetscReal *grtol, PetscReal *gttol)
1371: {
1372:   PetscFunctionBegin;
1374:   if (gatol) *gatol = tao->gatol;
1375:   if (grtol) *grtol = tao->grtol;
1376:   if (gttol) *gttol = tao->gttol;
1377:   PetscFunctionReturn(PETSC_SUCCESS);
1378: }

1380: /*@
1381:   TaoGetKSP - Gets the linear solver used by the optimization solver.

1383:   Not Collective

1385:   Input Parameter:
1386: . tao - the `Tao` solver

1388:   Output Parameter:
1389: . ksp - the `KSP` linear solver used in the optimization solver

1391:   Level: intermediate

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

1404: /*@
1405:   TaoGetLinearSolveIterations - Gets the total number of linear iterations
1406:   used by the `Tao` solver

1408:   Not Collective

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

1413:   Output Parameter:
1414: . lits - number of linear iterations

1416:   Level: intermediate

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

1421: .seealso: [](ch_tao), `Tao`, `TaoGetKSP()`
1422: @*/
1423: PetscErrorCode TaoGetLinearSolveIterations(Tao tao, PetscInt *lits)
1424: {
1425:   PetscFunctionBegin;
1427:   PetscAssertPointer(lits, 2);
1428:   *lits = tao->ksp_tot_its;
1429:   PetscFunctionReturn(PETSC_SUCCESS);
1430: }

1432: /*@
1433:   TaoGetLineSearch - Gets the line search used by the optimization solver.

1435:   Not Collective

1437:   Input Parameter:
1438: . tao - the `Tao` solver

1440:   Output Parameter:
1441: . ls - the line search used in the optimization solver

1443:   Level: intermediate

1445: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`
1446: @*/
1447: PetscErrorCode TaoGetLineSearch(Tao tao, TaoLineSearch *ls)
1448: {
1449:   PetscFunctionBegin;
1451:   PetscAssertPointer(ls, 2);
1452:   *ls = tao->linesearch;
1453:   PetscFunctionReturn(PETSC_SUCCESS);
1454: }

1456: /*@
1457:   TaoAddLineSearchCounts - Adds the number of function evaluations spent
1458:   in the line search to the running total.

1460:   Input Parameters:
1461: . tao - the `Tao` solver

1463:   Level: developer

1465: .seealso: [](ch_tao), `Tao`, `TaoGetLineSearch()`, `TaoLineSearchApply()`
1466: @*/
1467: PetscErrorCode TaoAddLineSearchCounts(Tao tao)
1468: {
1469:   PetscBool flg;
1470:   PetscInt  nfeval, ngeval, nfgeval;

1472:   PetscFunctionBegin;
1474:   if (tao->linesearch) {
1475:     PetscCall(TaoLineSearchIsUsingTaoRoutines(tao->linesearch, &flg));
1476:     if (!flg) {
1477:       PetscCall(TaoLineSearchGetNumberFunctionEvaluations(tao->linesearch, &nfeval, &ngeval, &nfgeval));
1478:       tao->objective_term.term->nobj += nfeval;
1479:       tao->objective_term.term->ngrad += ngeval;
1480:       tao->objective_term.term->nobjgrad += nfgeval;
1481:     }
1482:   }
1483:   PetscFunctionReturn(PETSC_SUCCESS);
1484: }

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

1489:   Not Collective

1491:   Input Parameter:
1492: . tao - the `Tao` context

1494:   Output Parameter:
1495: . X - the current solution

1497:   Level: intermediate

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

1502: .seealso: [](ch_tao), `Tao`, `TaoSetSolution()`, `TaoSolve()`
1503: @*/
1504: PetscErrorCode TaoGetSolution(Tao tao, Vec *X)
1505: {
1506:   PetscFunctionBegin;
1508:   PetscAssertPointer(X, 2);
1509:   *X = tao->solution;
1510:   PetscFunctionReturn(PETSC_SUCCESS);
1511: }

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

1518:   Collective

1520:   Input Parameter:
1521: . tao - the `Tao` context

1523:   Level: developer

1525:   Note:
1526:   This function does not reset the statistics of internal `TaoTerm`

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

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

1554:   Logically Collective

1556:   Input Parameters:
1557: + tao  - The `Tao` solver
1558: . func - The function
1559: - ctx  - The update function context

1561:   Calling sequence of `func`:
1562: + tao - The optimizer context
1563: . it  - The current iteration index
1564: - ctx - The update context

1566:   Level: advanced

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

1572: .seealso: [](ch_tao), `Tao`, `TaoSolve()`
1573: @*/
1574: PetscErrorCode TaoSetUpdate(Tao tao, PetscErrorCode (*func)(Tao tao, PetscInt it, PetscCtx ctx), PetscCtx ctx)
1575: {
1576:   PetscFunctionBegin;
1578:   tao->ops->update = func;
1579:   tao->user_update = ctx;
1580:   PetscFunctionReturn(PETSC_SUCCESS);
1581: }

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

1588:   Logically Collective

1590:   Input Parameters:
1591: + tao  - the `Tao` object
1592: . conv - the routine to test for convergence
1593: - ctx  - [optional] context for private data for the convergence routine (may be `NULL`)

1595:   Calling sequence of `conv`:
1596: + tao - the `Tao` object
1597: - ctx - [optional] convergence context

1599:   Level: advanced

1601:   Note:
1602:   The new convergence testing routine should call `TaoSetConvergedReason()`.

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

1615: /*@C
1616:   TaoMonitorSet - Sets an additional function that is to be used at every
1617:   iteration of the solver to display the iteration's
1618:   progress.

1620:   Logically Collective

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

1628:   Calling sequence of `func`:
1629: + tao - the `Tao` solver context
1630: - ctx - [optional] monitoring context

1632:   Level: intermediate

1634:   Notes:
1635:   See `TaoSetFromOptions()` for a monitoring options.

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

1641:   Fortran Notes:
1642:   Only one monitor function may be set

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

1654:     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))(PetscVoidFn *)func, ctx, dest, (PetscErrorCode (*)(void))(PetscVoidFn *)tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i], &identical));
1655:     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
1656:   }
1657:   tao->monitor[tao->numbermonitors]        = func;
1658:   tao->monitorcontext[tao->numbermonitors] = ctx;
1659:   tao->monitordestroy[tao->numbermonitors] = dest;
1660:   ++tao->numbermonitors;
1661:   PetscFunctionReturn(PETSC_SUCCESS);
1662: }

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

1667:   Logically Collective

1669:   Input Parameter:
1670: . tao - the `Tao` solver context

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

1677:   Level: advanced

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

1682: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1683: @*/
1684: PetscErrorCode TaoMonitorCancel(Tao tao)
1685: {
1686:   PetscFunctionBegin;
1688:   for (PetscInt i = 0; i < tao->numbermonitors; i++) {
1689:     if (tao->monitordestroy[i]) PetscCall((*tao->monitordestroy[i])(&tao->monitorcontext[i]));
1690:   }
1691:   tao->numbermonitors = 0;
1692:   PetscFunctionReturn(PETSC_SUCCESS);
1693: }

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

1698:   Collective

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

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

1707:   Level: advanced

1709:   Note:
1710:   This monitor prints the function value and gradient
1711:   norm at each iteration.

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

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

1725:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1726:   PetscCall(PetscViewerPushFormat(viewer, vf->format));
1727:   if (isascii) {
1728:     PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));

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

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

1751:   Collective

1753:   Input Parameters:
1754: + tao - the `Tao` context
1755: - vf  - `PetscViewerAndFormat` context

1757:   Options Database Key:
1758: . -tao_monitor_globalization - turn on monitoring with globalization information

1760:   Level: advanced

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

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

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

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

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

1805:   Collective

1807:   Input Parameters:
1808: + tao - the `Tao` context
1809: - vf  - `PetscViewerAndFormat` context

1811:   Options Database Key:
1812: . -tao_monitor_short - turn on default short monitoring

1814:   Level: advanced

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

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

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

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

1859: /*@
1860:   TaoMonitorConstraintNorm - same as `TaoMonitorDefault()` except
1861:   it prints the norm of the constraint function.

1863:   Collective

1865:   Input Parameters:
1866: + tao - the `Tao` context
1867: - vf  - `PetscViewerAndFormat` context

1869:   Options Database Key:
1870: . -tao_monitor_constraint_norm - monitor the constraints

1872:   Level: advanced

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

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

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

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

1904:   Collective

1906:   Input Parameters:
1907: + tao - the `Tao` context
1908: - vf  - `PetscViewerAndFormat` context

1910:   Options Database Key:
1911: . -tao_monitor_solution - view the solution

1913:   Level: advanced

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

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

1931:   Collective

1933:   Input Parameters:
1934: + tao - the `Tao` context
1935: - vf  - `PetscViewerAndFormat` context

1937:   Options Database Key:
1938: . -tao_monitor_gradient - view the gradient at each iteration

1940:   Level: advanced

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

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

1958:   Collective

1960:   Input Parameters:
1961: + tao - the `Tao` context
1962: - vf  - `PetscViewerAndFormat` context

1964:   Options Database Key:
1965: . -tao_monitor_step - view the step vector at each iteration

1967:   Level: advanced

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

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

1985:   Collective

1987:   Input Parameters:
1988: + tao - the `Tao` context
1989: - ctx - `TaoMonitorDraw` context

1991:   Options Database Key:
1992: . -tao_monitor_solution_draw - draw the solution at each iteration

1994:   Level: advanced

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

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

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

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

2017:   Collective

2019:   Input Parameters:
2020: + tao - the `Tao` context
2021: - ctx - `PetscViewer` context

2023:   Options Database Key:
2024: . -tao_monitor_gradient_draw - draw the gradient at each iteration

2026:   Level: advanced

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

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

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

2044:   Collective

2046:   Input Parameters:
2047: + tao - the `Tao` context
2048: - ctx - the `PetscViewer` context

2050:   Options Database Key:
2051: . -tao_monitor_step_draw - draw the step direction at each iteration

2053:   Level: advanced

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

2061:   PetscFunctionBegin;
2064:   PetscCall(VecView(tao->stepdirection, viewer));
2065:   PetscFunctionReturn(PETSC_SUCCESS);
2066: }

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

2071:   Collective

2073:   Input Parameters:
2074: + tao - the `Tao` context
2075: - vf  - `PetscViewerAndFormat` context

2077:   Options Database Key:
2078: . -tao_monitor_ls_residual - view the residual at each iteration

2080:   Level: advanced

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

2095: /*@
2096:   TaoDefaultConvergenceTest - Determines whether the solver should continue iterating
2097:   or terminate.

2099:   Collective

2101:   Input Parameters:
2102: + tao   - the `Tao` context
2103: - dummy - unused dummy context

2105:   Level: developer

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

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

2126:   PetscFunctionBegin;
2128:   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);

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

2165: /*@
2166:   TaoSetOptionsPrefix - Sets the prefix used for searching for all
2167:   Tao options in the database.

2169:   Logically Collective

2171:   Input Parameters:
2172: + tao - the `Tao` context
2173: - p   - the prefix string to prepend to all Tao option requests

2175:   Level: advanced

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

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

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

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

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

2213:   Logically Collective

2215:   Input Parameters:
2216: + tao - the `Tao` solver context
2217: - p   - the prefix string to prepend to all `Tao` option requests

2219:   Level: advanced

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

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

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

2244: /*@
2245:   TaoGetOptionsPrefix - Gets the prefix used for searching for all
2246:   Tao options in the database

2248:   Not Collective

2250:   Input Parameter:
2251: . tao - the `Tao` context

2253:   Output Parameter:
2254: . p - pointer to the prefix string used is returned

2256:   Level: advanced

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

2268: /*@
2269:   TaoSetType - Sets the `TaoType` for the minimization solver.

2271:   Collective

2273:   Input Parameters:
2274: + tao  - the `Tao` solver context
2275: - type - a known method

2277:   Options Database Key:
2278: . -tao_type type - Sets the method; see `TaoType`

2280:   Level: intermediate

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

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

2294:   PetscFunctionBegin;

2297:   PetscCall(PetscObjectTypeCompare((PetscObject)tao, type, &issame));
2298:   if (issame) PetscFunctionReturn(PETSC_SUCCESS);

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

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

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

2317:   tao->setupcalled           = PETSC_FALSE;
2318:   tao->uses_gradient         = PETSC_FALSE;
2319:   tao->uses_hessian_matrices = PETSC_FALSE;

2321:   PetscCall((*create_xxx)(tao));
2322:   PetscCall(PetscObjectChangeTypeName((PetscObject)tao, type));
2323:   PetscFunctionReturn(PETSC_SUCCESS);
2324: }

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

2329:   Not Collective, No Fortran Support

2331:   Input Parameters:
2332: + sname - name of a new user-defined solver
2333: - func  - routine to create `TaoType` specific method context

2335:   Calling sequence of `func`:
2336: . tao - the `Tao` object to be created

2338:   Example Usage:
2339: .vb
2340:    TaoRegister("my_solver", MySolverCreate);
2341: .ve

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

2352:   Level: advanced

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

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

2367: /*@C
2368:   TaoRegisterDestroy - Frees the list of minimization solvers that were
2369:   registered by `TaoRegister()`.

2371:   Not Collective

2373:   Level: advanced

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

2385: /*@
2386:   TaoGetIterationNumber - Gets the number of `TaoSolve()` iterations completed
2387:   at this time.

2389:   Not Collective

2391:   Input Parameter:
2392: . tao - the `Tao` context

2394:   Output Parameter:
2395: . iter - iteration number

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

2400:   Level: intermediate

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

2413: /*@
2414:   TaoGetResidualNorm - Gets the current value of the norm of the residual (gradient)
2415:   at this time.

2417:   Not Collective

2419:   Input Parameter:
2420: . tao - the `Tao` context

2422:   Output Parameter:
2423: . value - the current value

2425:   Level: intermediate

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

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

2442: /*@
2443:   TaoSetIterationNumber - Sets the current iteration number.

2445:   Logically Collective

2447:   Input Parameters:
2448: + tao  - the `Tao` context
2449: - iter - iteration number

2451:   Level: developer

2453: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
2454: @*/
2455: PetscErrorCode TaoSetIterationNumber(Tao tao, PetscInt iter)
2456: {
2457:   PetscFunctionBegin;
2460:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2461:   tao->niter = iter;
2462:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2463:   PetscFunctionReturn(PETSC_SUCCESS);
2464: }

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

2471:   Not Collective

2473:   Input Parameter:
2474: . tao - the `Tao` context

2476:   Output Parameter:
2477: . iter - number of iterations

2479:   Level: intermediate

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

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

2496: /*@
2497:   TaoSetTotalIterationNumber - Sets the current total iteration number.

2499:   Logically Collective

2501:   Input Parameters:
2502: + tao  - the `Tao` context
2503: - iter - the iteration number

2505:   Level: developer

2507: .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
2508: @*/
2509: PetscErrorCode TaoSetTotalIterationNumber(Tao tao, PetscInt iter)
2510: {
2511:   PetscFunctionBegin;
2514:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2515:   tao->ntotalits = iter;
2516:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2517:   PetscFunctionReturn(PETSC_SUCCESS);
2518: }

2520: /*@
2521:   TaoSetConvergedReason - Sets the termination flag on a `Tao` object

2523:   Logically Collective

2525:   Input Parameters:
2526: + tao    - the `Tao` context
2527: - reason - the `TaoConvergedReason`

2529:   Level: intermediate

2531: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`
2532: @*/
2533: PetscErrorCode TaoSetConvergedReason(Tao tao, TaoConvergedReason reason)
2534: {
2535:   PetscFunctionBegin;
2538:   tao->reason = reason;
2539:   PetscFunctionReturn(PETSC_SUCCESS);
2540: }

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

2545:   Not Collective

2547:   Input Parameter:
2548: . tao - the `Tao` solver context

2550:   Output Parameter:
2551: . reason - value of `TaoConvergedReason`

2553:   Level: intermediate

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

2566: /*@
2567:   TaoGetSolutionStatus - Get the current iterate, objective value,
2568:   residual, infeasibility, and termination from a `Tao` object

2570:   Not Collective

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

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

2583:   Level: intermediate

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

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

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

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

2608:   Not Collective

2610:   Input Parameter:
2611: . tao - the `Tao` solver context

2613:   Output Parameter:
2614: . type - the `TaoType`

2616:   Level: intermediate

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

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

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

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

2646:   Level: developer

2648: .seealso: [](ch_tao), `Tao`, `TaoGetConvergedReason()`, `TaoMonitorDefault()`, `TaoMonitorSet()`
2649: @*/
2650: PetscErrorCode TaoMonitor(Tao tao, PetscInt its, PetscReal f, PetscReal res, PetscReal cnorm, PetscReal steplength)
2651: {
2652:   PetscFunctionBegin;
2654:   tao->fc       = f;
2655:   tao->residual = res;
2656:   tao->cnorm    = cnorm;
2657:   tao->step     = steplength;
2658:   if (!its) {
2659:     tao->cnorm0 = cnorm;
2660:     tao->gnorm0 = res;
2661:   }
2662:   PetscCall(VecLockReadPush(tao->solution));
2663:   for (PetscInt i = 0; i < tao->numbermonitors; i++) PetscCall((*tao->monitor[i])(tao, tao->monitorcontext[i]));
2664:   PetscCall(VecLockReadPop(tao->solution));
2665:   PetscFunctionReturn(PETSC_SUCCESS);
2666: }

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

2671:   Logically Collective

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

2683:   Level: intermediate

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

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

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

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

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

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

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

2728:   Collective

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

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

2740:   Level: advanced

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

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

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

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

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

2775:   Logically Collective

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

2781:   Level: intermediate

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

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

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

2801:   Not Collective

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

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

2809:   Level: intermediate

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

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

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

2831:   Collective

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

2837:   Level: beginner

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

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

2857:   Not Collective

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

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

2865:   Level: beginner

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

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

2881:   Collective

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

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

2891:   Level: advanced

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

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

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

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

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

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

2927:   Collective

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

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

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

2946:   Level: intermediate

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

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

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

2967:   Collective

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

2972:   Level: intermediate

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

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

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

2993:   Not collective

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

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

3004:   Level: intermediate

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

3011:   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.

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

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

3032:   Collective

3034:   Input Parameters:
3035: + tao    - a `Tao` solver context
3036: . 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.)
3037: . scale  - scaling coefficient for the new term
3038: . term   - the real-valued function defining the new term
3039: . params - (optional) parameters for the new term.  It is up to each implementation of `TaoTerm` to determine how it behaves when parameters are omitted.
3040: - map    - (optional) a map from the `tao` solution space to the `term` solution space; if `NULL` the map is assumed to be the identity

3042:   Level: beginner

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

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

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

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

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

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

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

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

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

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

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

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

3175: /*@
3176:   TaoSetDM - Sets the `DM` that may be used by some `TAO` solvers or their underlying solvers and preconditioners

3178:   Logically Collective

3180:   Input Parameters:
3181: + tao - the nonlinear solver context
3182: - dm  - the `DM`, cannot be `NULL`

3184:   Level: intermediate

3186:   Note:
3187:   A `DM` can only be used for solving one problem at a time because information about the problem is stored on the `DM`,
3188:   even when not using interfaces like `DMSNESSetFunction()`.  Use `DMClone()` to get a distinct `DM` when solving different
3189:   problems using the same function space.

3191: .seealso: [](ch_snes), `DM`, `TAO`, `TaoGetDM()`, `SNESSetDM()`, `SNESGetDM()`, `KSPSetDM()`, `KSPGetDM()`
3192: @*/
3193: PetscErrorCode TaoSetDM(Tao tao, DM dm)
3194: {
3195:   KSP ksp;

3197:   PetscFunctionBegin;
3200:   PetscCall(PetscObjectReference((PetscObject)dm));
3201:   PetscCall(DMDestroy(&tao->dm));
3202:   tao->dm = dm;

3204:   PetscCall(TaoGetKSP(tao, &ksp));
3205:   if (ksp) {
3206:     PetscCall(KSPSetDM(ksp, dm));
3207:     PetscCall(KSPSetDMActive(ksp, KSP_DMACTIVE_ALL, PETSC_FALSE));
3208:   }
3209:   PetscFunctionReturn(PETSC_SUCCESS);
3210: }

3212: /*@
3213:   TaoGetDM - Gets the `DM` that may be used by some `TAO` solvers or their underlying solvers and preconditioners

3215:   Not Collective but `dm` obtained is parallel on `tao`

3217:   Input Parameter:
3218: . tao - the `TAO` context

3220:   Output Parameter:
3221: . dm - the `DM`

3223:   Level: intermediate

3225: .seealso: [](ch_snes), `DM`, `TAO`, `TaoSetDM()`, `SNESSetDM()`, `SNESGetDM()`, `KSPSetDM()`, `KSPGetDM()`
3226: @*/
3227: PetscErrorCode TaoGetDM(Tao tao, DM *dm)
3228: {
3229:   PetscFunctionBegin;
3231:   PetscAssertPointer(dm, 2);
3232:   if (!tao->dm) PetscCall(DMShellCreate(PetscObjectComm((PetscObject)tao), &tao->dm));
3233:   *dm = tao->dm;
3234:   PetscFunctionReturn(PETSC_SUCCESS);
3235: }