Actual source code: taosolver.c
1: #include <petsc/private/taoimpl.h>
2: #include <petsc/private/snesimpl.h>
4: PetscBool TaoRegisterAllCalled = PETSC_FALSE;
5: PetscFunctionList TaoList = NULL;
7: PetscClassId TAO_CLASSID = 0;
9: PetscLogEvent TAO_Solve;
10: PetscLogEvent TAO_ResidualEval;
11: PetscLogEvent TAO_JacobianEval;
12: PetscLogEvent TAO_ConstraintsEval;
14: const char *TaoSubSetTypes[] = {"subvec", "mask", "matrixfree", "TaoSubSetType", "TAO_SUBSET_", NULL};
16: struct _n_TaoMonitorDrawCtx {
17: PetscViewer viewer;
18: PetscInt howoften; /* when > 0 uses iteration % howoften, when negative only final solution plotted */
19: };
21: static PetscErrorCode KSPPreSolve_TAOEW_Private(KSP ksp, Vec b, Vec x, PetscCtx ctx)
22: {
23: Tao tao = (Tao)ctx;
24: SNES snes_ewdummy = tao->snes_ewdummy;
26: PetscFunctionBegin;
27: if (!snes_ewdummy) PetscFunctionReturn(PETSC_SUCCESS);
28: /* populate snes_ewdummy struct values used in KSPPreSolve_SNESEW */
29: snes_ewdummy->vec_func = b;
30: snes_ewdummy->rtol = tao->gttol;
31: snes_ewdummy->iter = tao->niter;
32: PetscCall(VecNorm(b, NORM_2, &snes_ewdummy->norm));
33: PetscCall(KSPPreSolve_SNESEW(ksp, b, x, snes_ewdummy));
34: snes_ewdummy->vec_func = NULL;
35: PetscFunctionReturn(PETSC_SUCCESS);
36: }
38: static PetscErrorCode KSPPostSolve_TAOEW_Private(KSP ksp, Vec b, Vec x, PetscCtx ctx)
39: {
40: Tao tao = (Tao)ctx;
41: SNES snes_ewdummy = tao->snes_ewdummy;
43: PetscFunctionBegin;
44: if (!snes_ewdummy) PetscFunctionReturn(PETSC_SUCCESS);
45: PetscCall(KSPPostSolve_SNESEW(ksp, b, x, snes_ewdummy));
46: PetscFunctionReturn(PETSC_SUCCESS);
47: }
49: static PetscErrorCode TaoSetUpEW_Private(Tao tao)
50: {
51: SNESKSPEW *kctx;
52: const char *ewprefix;
54: PetscFunctionBegin;
55: if (!tao->ksp) PetscFunctionReturn(PETSC_SUCCESS);
56: if (tao->ksp_ewconv) {
57: if (!tao->snes_ewdummy) PetscCall(SNESCreate(PetscObjectComm((PetscObject)tao), &tao->snes_ewdummy));
58: tao->snes_ewdummy->ksp_ewconv = PETSC_TRUE;
59: PetscCall(KSPSetPreSolve(tao->ksp, KSPPreSolve_TAOEW_Private, tao));
60: PetscCall(KSPSetPostSolve(tao->ksp, KSPPostSolve_TAOEW_Private, tao));
62: PetscCall(KSPGetOptionsPrefix(tao->ksp, &ewprefix));
63: kctx = (SNESKSPEW *)tao->snes_ewdummy->kspconvctx;
64: PetscCall(SNESEWSetFromOptions_Private(kctx, PETSC_FALSE, PetscObjectComm((PetscObject)tao), ewprefix));
65: } else PetscCall(SNESDestroy(&tao->snes_ewdummy));
66: PetscFunctionReturn(PETSC_SUCCESS);
67: }
69: /*@
70: TaoParametersInitialize - Sets all the parameters in `tao` to their default value (when `TaoCreate()` was called) if they
71: currently contain default values. Default values are the parameter values when the object's type is set.
73: Collective
75: Input Parameter:
76: . tao - the `Tao` object
78: Level: developer
80: Developer Note:
81: This is called by all the `TaoCreate_XXX()` routines.
83: .seealso: [](ch_snes), `Tao`, `TaoSolve()`, `TaoDestroy()`,
84: `PetscObjectParameterSetDefault()`
85: @*/
86: PetscErrorCode TaoParametersInitialize(Tao tao)
87: {
88: PetscObjectParameterSetDefault(tao, max_it, 10000);
89: PetscObjectParameterSetDefault(tao, max_funcs, PETSC_UNLIMITED);
90: PetscObjectParameterSetDefault(tao, gatol, PetscDefined(USE_REAL_SINGLE) ? 1e-5 : 1e-8);
91: PetscObjectParameterSetDefault(tao, grtol, PetscDefined(USE_REAL_SINGLE) ? 1e-5 : 1e-8);
92: PetscObjectParameterSetDefault(tao, crtol, PetscDefined(USE_REAL_SINGLE) ? 1e-5 : 1e-8);
93: PetscObjectParameterSetDefault(tao, catol, PetscDefined(USE_REAL_SINGLE) ? 1e-5 : 1e-8);
94: PetscObjectParameterSetDefault(tao, gttol, 0.0);
95: PetscObjectParameterSetDefault(tao, steptol, 0.0);
96: PetscObjectParameterSetDefault(tao, fmin, PETSC_NINFINITY);
97: PetscObjectParameterSetDefault(tao, trust0, PETSC_INFINITY);
98: return PETSC_SUCCESS;
99: }
101: /*@
102: TaoCreate - Creates a Tao solver
104: Collective
106: Input Parameter:
107: . comm - MPI communicator
109: Output Parameter:
110: . newtao - the new `Tao` context
112: Options Database Key:
113: . -tao_type - select which method Tao should use
115: Level: beginner
117: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoDestroy()`, `TaoSetFromOptions()`, `TaoSetType()`
118: @*/
119: PetscErrorCode TaoCreate(MPI_Comm comm, Tao *newtao)
120: {
121: Tao tao;
123: PetscFunctionBegin;
124: PetscAssertPointer(newtao, 2);
125: PetscCall(TaoInitializePackage());
126: PetscCall(TaoLineSearchInitializePackage());
128: PetscCall(PetscHeaderCreate(tao, TAO_CLASSID, "Tao", "Optimization solver", "Tao", comm, TaoDestroy, TaoView));
129: tao->ops->convergencetest = TaoDefaultConvergenceTest;
131: tao->hist_reset = PETSC_TRUE;
132: tao->term_set = PETSC_FALSE;
134: PetscCall(TaoTermCreateCallbacks(tao, &tao->callbacks));
135: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tao->callbacks, "callbacks_"));
136: PetscCall(TaoTermMappingSetData(&tao->objective_term, NULL, 1.0, tao->callbacks, NULL));
137: PetscCall(TaoResetStatistics(tao));
138: *newtao = tao;
139: PetscFunctionReturn(PETSC_SUCCESS);
140: }
142: /*@
143: TaoSolve - Solves an optimization problem min F(x) s.t. l <= x <= u
145: Collective
147: Input Parameter:
148: . tao - the `Tao` context
150: Level: beginner
152: Notes:
153: The user must set up the `Tao` object with calls to `TaoSetSolution()`, `TaoSetObjective()`, `TaoSetGradient()`, and (if using 2nd order method) `TaoSetHessian()`.
155: You should call `TaoGetConvergedReason()` or run with `-tao_converged_reason` to determine if the optimization algorithm actually succeeded or
156: why it failed.
158: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSetObjective()`, `TaoSetGradient()`, `TaoSetHessian()`, `TaoGetConvergedReason()`, `TaoSetUp()`
159: @*/
160: PetscErrorCode TaoSolve(Tao tao)
161: {
162: static PetscBool set = PETSC_FALSE;
164: PetscFunctionBegin;
166: PetscCall(PetscCitationsRegister("@TechReport{tao-user-ref,\n"
167: "title = {Toolkit for Advanced Optimization (TAO) Users Manual},\n"
168: "author = {Todd Munson and Jason Sarich and Stefan Wild and Steve Benson and Lois Curfman McInnes},\n"
169: "Institution = {Argonne National Laboratory},\n"
170: "Year = 2014,\n"
171: "Number = {ANL/MCS-TM-322 - Revision 3.5},\n"
172: "url = {https://www.mcs.anl.gov/research/projects/tao/}\n}\n",
173: &set));
174: tao->header_printed = PETSC_FALSE;
175: PetscCall(TaoSetUp(tao));
176: PetscCall(TaoResetStatistics(tao));
177: if (tao->linesearch) PetscCall(TaoLineSearchReset(tao->linesearch));
179: PetscCall(PetscLogEventBegin(TAO_Solve, tao, 0, 0, 0));
180: PetscTryTypeMethod(tao, solve);
181: PetscCall(PetscLogEventEnd(TAO_Solve, tao, 0, 0, 0));
183: PetscCall(VecViewFromOptions(tao->solution, (PetscObject)tao, "-tao_view_solution"));
185: tao->ntotalits += tao->niter;
187: if (tao->printreason) {
188: PetscViewer viewer = PETSC_VIEWER_STDOUT_(((PetscObject)tao)->comm);
190: PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)tao)->tablevel));
191: if (tao->reason > 0) {
192: if (((PetscObject)tao)->prefix) {
193: PetscCall(PetscViewerASCIIPrintf(viewer, "TAO %s solve converged due to %s iterations %" PetscInt_FMT "\n", ((PetscObject)tao)->prefix, TaoConvergedReasons[tao->reason], tao->niter));
194: } else {
195: PetscCall(PetscViewerASCIIPrintf(viewer, "TAO solve converged due to %s iterations %" PetscInt_FMT "\n", TaoConvergedReasons[tao->reason], tao->niter));
196: }
197: } else {
198: if (((PetscObject)tao)->prefix) {
199: PetscCall(PetscViewerASCIIPrintf(viewer, "TAO %s solve did not converge due to %s iteration %" PetscInt_FMT "\n", ((PetscObject)tao)->prefix, TaoConvergedReasons[tao->reason], tao->niter));
200: } else {
201: PetscCall(PetscViewerASCIIPrintf(viewer, "TAO solve did not converge due to %s iteration %" PetscInt_FMT "\n", TaoConvergedReasons[tao->reason], tao->niter));
202: }
203: }
204: PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)tao)->tablevel));
205: }
206: PetscCall(TaoViewFromOptions(tao, NULL, "-tao_view"));
207: PetscFunctionReturn(PETSC_SUCCESS);
208: }
210: /*@
211: TaoSetUp - Sets up the internal data structures for the later use
212: of a Tao solver
214: Collective
216: Input Parameter:
217: . tao - the `Tao` context
219: Level: advanced
221: Note:
222: The user will not need to explicitly call `TaoSetUp()`, as it will
223: automatically be called in `TaoSolve()`. However, if the user
224: desires to call it explicitly, it should come after `TaoCreate()`
225: and any TaoSetSomething() routines, but before `TaoSolve()`.
227: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
228: @*/
229: PetscErrorCode TaoSetUp(Tao tao)
230: {
231: PetscFunctionBegin;
233: if (tao->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
234: PetscCall(TaoSetUpEW_Private(tao));
235: PetscCall(TaoTermMappingSetUp(&tao->objective_term));
236: if (!tao->solution) PetscCall(TaoTermMappingCreateSolutionVec(&tao->objective_term, &tao->solution));
237: PetscCheck(tao->solution, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoSetSolution()");
238: if (tao->uses_gradient && !tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
239: if (tao->uses_hessian_matrices) {
240: // TaoSetHessian has been called, but as terms have been added,
241: // subterms' Hessian and PtAP routines, if needed, have to be created
242: // TODO Function to set TAOTERMSUM's Hessian.
243: if (!tao->hessian) {
244: PetscBool is_defined;
246: // TAOTERMSUM's Hessian will follow layout and type of first term's Hessian
247: PetscCall(TaoTermIsCreateHessianMatricesDefined(tao->objective_term.term, &is_defined));
248: if (is_defined) PetscCall(TaoTermMappingCreateHessianMatrices(&tao->objective_term, &tao->hessian, &tao->hessian_pre));
249: }
250: PetscCheck(tao->hessian, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoSetHessian()");
251: }
252: PetscTryTypeMethod(tao, setup);
253: tao->setupcalled = PETSC_TRUE;
254: PetscFunctionReturn(PETSC_SUCCESS);
255: }
257: /*@
258: TaoDestroy - Destroys the `Tao` context that was created with `TaoCreate()`
260: Collective
262: Input Parameter:
263: . tao - the `Tao` context
265: Level: beginner
267: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
268: @*/
269: PetscErrorCode TaoDestroy(Tao *tao)
270: {
271: PetscFunctionBegin;
272: if (!*tao) PetscFunctionReturn(PETSC_SUCCESS);
274: if (--((PetscObject)*tao)->refct > 0) {
275: *tao = NULL;
276: PetscFunctionReturn(PETSC_SUCCESS);
277: }
279: PetscTryTypeMethod(*tao, destroy);
280: PetscCall(TaoTermMappingReset(&(*tao)->objective_term));
281: PetscCall(VecDestroy(&(*tao)->objective_parameters));
282: PetscCall(TaoTermDestroy(&(*tao)->callbacks));
283: PetscCall(KSPDestroy(&(*tao)->ksp));
284: PetscCall(SNESDestroy(&(*tao)->snes_ewdummy));
285: PetscCall(TaoLineSearchDestroy(&(*tao)->linesearch));
287: if ((*tao)->ops->convergencedestroy) {
288: PetscCall((*(*tao)->ops->convergencedestroy)((*tao)->cnvP));
289: PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
290: }
291: PetscCall(VecDestroy(&(*tao)->solution));
292: PetscCall(VecDestroy(&(*tao)->gradient));
293: PetscCall(VecDestroy(&(*tao)->ls_res));
295: if ((*tao)->gradient_norm) {
296: PetscCall(PetscObjectDereference((PetscObject)(*tao)->gradient_norm));
297: PetscCall(VecDestroy(&(*tao)->gradient_norm_tmp));
298: }
300: PetscCall(VecDestroy(&(*tao)->XL));
301: PetscCall(VecDestroy(&(*tao)->XU));
302: PetscCall(VecDestroy(&(*tao)->IL));
303: PetscCall(VecDestroy(&(*tao)->IU));
304: PetscCall(VecDestroy(&(*tao)->DE));
305: PetscCall(VecDestroy(&(*tao)->DI));
306: PetscCall(VecDestroy(&(*tao)->constraints));
307: PetscCall(VecDestroy(&(*tao)->constraints_equality));
308: PetscCall(VecDestroy(&(*tao)->constraints_inequality));
309: PetscCall(VecDestroy(&(*tao)->stepdirection));
310: PetscCall(MatDestroy(&(*tao)->hessian_pre));
311: PetscCall(MatDestroy(&(*tao)->hessian));
312: PetscCall(MatDestroy(&(*tao)->ls_jac));
313: PetscCall(MatDestroy(&(*tao)->ls_jac_pre));
314: PetscCall(MatDestroy(&(*tao)->jacobian_pre));
315: PetscCall(MatDestroy(&(*tao)->jacobian));
316: PetscCall(MatDestroy(&(*tao)->jacobian_state_pre));
317: PetscCall(MatDestroy(&(*tao)->jacobian_state));
318: PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
319: PetscCall(MatDestroy(&(*tao)->jacobian_design));
320: PetscCall(MatDestroy(&(*tao)->jacobian_equality));
321: PetscCall(MatDestroy(&(*tao)->jacobian_equality_pre));
322: PetscCall(MatDestroy(&(*tao)->jacobian_inequality));
323: PetscCall(MatDestroy(&(*tao)->jacobian_inequality_pre));
324: PetscCall(ISDestroy(&(*tao)->state_is));
325: PetscCall(ISDestroy(&(*tao)->design_is));
326: PetscCall(VecDestroy(&(*tao)->res_weights_v));
327: PetscCall(TaoMonitorCancel(*tao));
328: if ((*tao)->hist_malloc) PetscCall(PetscFree4((*tao)->hist_obj, (*tao)->hist_resid, (*tao)->hist_cnorm, (*tao)->hist_lits));
329: if ((*tao)->res_weights_n) {
330: PetscCall(PetscFree((*tao)->res_weights_rows));
331: PetscCall(PetscFree((*tao)->res_weights_cols));
332: PetscCall(PetscFree((*tao)->res_weights_w));
333: }
334: PetscCall(PetscHeaderDestroy(tao));
335: PetscFunctionReturn(PETSC_SUCCESS);
336: }
338: /*@
339: TaoKSPSetUseEW - Sets `SNES` to use Eisenstat-Walker method {cite}`ew96` for computing relative tolerance for linear solvers.
341: Logically Collective
343: Input Parameters:
344: + tao - Tao context
345: - flag - `PETSC_TRUE` or `PETSC_FALSE`
347: Level: advanced
349: Note:
350: See `SNESKSPSetUseEW()` for customization details.
352: .seealso: [](ch_tao), `Tao`, `SNESKSPSetUseEW()`
353: @*/
354: PetscErrorCode TaoKSPSetUseEW(Tao tao, PetscBool flag)
355: {
356: PetscFunctionBegin;
359: tao->ksp_ewconv = flag;
360: PetscFunctionReturn(PETSC_SUCCESS);
361: }
363: /*@C
364: TaoMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
366: Collective
368: Input Parameters:
369: + tao - `Tao` object you wish to monitor
370: . name - the monitor type one is seeking
371: . help - message indicating what monitoring is done
372: . manual - manual page for the monitor
373: - monitor - the monitor function, this must use a `PetscViewerFormat` as its context
375: Level: developer
377: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `PetscOptionsCreateViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`,
378: `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()`,
379: `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`,
380: `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`,
381: `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`,
382: `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`,
383: `PetscOptionsFList()`, `PetscOptionsEList()`
384: @*/
385: PetscErrorCode TaoMonitorSetFromOptions(Tao tao, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(Tao, PetscViewerAndFormat *))
386: {
387: PetscViewer viewer;
388: PetscViewerFormat format;
389: PetscBool flg;
391: PetscFunctionBegin;
392: PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)tao), ((PetscObject)tao)->options, ((PetscObject)tao)->prefix, name, &viewer, &format, &flg));
393: if (flg) {
394: PetscViewerAndFormat *vf;
395: char interval_key[1024];
397: PetscCall(PetscSNPrintf(interval_key, sizeof interval_key, "%s_interval", name));
398: PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
399: vf->view_interval = 1;
400: PetscCall(PetscOptionsGetInt(((PetscObject)tao)->options, ((PetscObject)tao)->prefix, interval_key, &vf->view_interval, NULL));
402: PetscCall(PetscViewerDestroy(&viewer));
403: PetscCall(TaoMonitorSet(tao, (PetscErrorCode (*)(Tao, PetscCtx))monitor, vf, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy));
404: }
405: PetscFunctionReturn(PETSC_SUCCESS);
406: }
408: /*@
409: TaoSetFromOptions - Sets various Tao parameters from the options database
411: Collective
413: Input Parameter:
414: . tao - the `Tao` solver context
416: Options Database Keys:
417: + -tao_type type - The algorithm that Tao uses (lmvm, nls, etc.). See `TAOType`
418: . -tao_gatol gatol - absolute error tolerance for ||gradient||
419: . -tao_grtol grtol - relative error tolerance for ||gradient||
420: . -tao_gttol gttol - reduction of ||gradient|| relative to initial gradient
421: . -tao_max_it max - sets maximum number of iterations
422: . -tao_max_funcs max - sets maximum number of function evaluations
423: . -tao_fmin fmin - stop if function value reaches fmin
424: . -tao_steptol tol - stop if trust region radius less than `tol`
425: . -tao_trust0 radius - initial trust region radius
426: . -tao_view_solution - view the solution at the end of the optimization process
427: . -tao_monitor - prints function value and residual norm at each iteration
428: . -tao_monitor_short - same as `-tao_monitor`, but truncates very small values
429: . -tao_monitor_constraint_norm - prints objective value, gradient, and constraint norm at each iteration
430: . -tao_monitor_globalization - prints information about the globalization at each iteration
431: . -tao_monitor_solution - prints solution vector at each iteration
432: . -tao_monitor_ls_residual - prints least-squares residual vector at each iteration
433: . -tao_monitor_step - prints step vector at each iteration
434: . -tao_monitor_gradient - prints gradient vector at each iteration
435: . -tao_monitor_solution_draw - graphically view solution vector at each iteration
436: . -tao_monitor_step_draw - graphically view step vector at each iteration
437: . -tao_monitor_gradient_draw - graphically view gradient at each iteration
438: . -tao_monitor_cancel - cancels all monitors (except those set with command line)
439: . -tao_fd_gradient - use gradient computed with finite differences
440: . -tao_fd_hessian - use hessian computed with finite differences
441: . -tao_mf_hessian - use matrix-free Hessian computed with finite differences. No `TaoTerm` support
442: . -tao_view - prints information about the Tao after solving
443: . -tao_converged_reason - prints the reason Tao stopped iterating
444: - -tao_add_terms - takes a comma-separated list of up to 16 options prefixes, a `TaoTerm` will be created for each and added to the objective function
446: Level: beginner
448: Notes:
449: To see all options, run your program with the `-help` option or consult the
450: user's manual. Should be called after `TaoCreate()` but before `TaoSolve()`.
452: The `-tao_add_terms` option accepts at most 16 prefixes.
454: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
455: @*/
456: PetscErrorCode TaoSetFromOptions(Tao tao)
457: {
458: TaoType default_type = TAOLMVM;
459: char type[256];
460: PetscBool flg, found;
461: MPI_Comm comm;
462: PetscReal catol, crtol, gatol, grtol, gttol;
464: PetscFunctionBegin;
466: PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
468: if (((PetscObject)tao)->type_name) default_type = ((PetscObject)tao)->type_name;
470: PetscObjectOptionsBegin((PetscObject)tao);
471: /* Check for type from options */
472: PetscCall(PetscOptionsFList("-tao_type", "Tao Solver type", "TaoSetType", TaoList, default_type, type, 256, &flg));
473: if (flg) PetscCall(TaoSetType(tao, type));
474: else if (!((PetscObject)tao)->type_name) PetscCall(TaoSetType(tao, default_type));
476: /* Tao solvers do not set the prefix, set it here if not yet done
477: We do it after SetType since solver may have been changed */
478: if (tao->linesearch) {
479: const char *prefix;
480: PetscCall(TaoLineSearchGetOptionsPrefix(tao->linesearch, &prefix));
481: if (!prefix) PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, ((PetscObject)tao)->prefix));
482: }
484: catol = tao->catol;
485: crtol = tao->crtol;
486: PetscCall(PetscOptionsReal("-tao_catol", "Stop if constraints violations within", "TaoSetConstraintTolerances", tao->catol, &catol, NULL));
487: PetscCall(PetscOptionsReal("-tao_crtol", "Stop if relative constraint violations within", "TaoSetConstraintTolerances", tao->crtol, &crtol, NULL));
488: PetscCall(TaoSetConstraintTolerances(tao, catol, crtol));
490: gatol = tao->gatol;
491: grtol = tao->grtol;
492: gttol = tao->gttol;
493: PetscCall(PetscOptionsReal("-tao_gatol", "Stop if norm of gradient less than", "TaoSetTolerances", tao->gatol, &gatol, NULL));
494: PetscCall(PetscOptionsReal("-tao_grtol", "Stop if norm of gradient divided by the function value is less than", "TaoSetTolerances", tao->grtol, &grtol, NULL));
495: PetscCall(PetscOptionsReal("-tao_gttol", "Stop if the norm of the gradient is less than the norm of the initial gradient times tol", "TaoSetTolerances", tao->gttol, >tol, NULL));
496: PetscCall(TaoSetTolerances(tao, gatol, grtol, gttol));
498: PetscCall(PetscOptionsInt("-tao_max_it", "Stop if iteration number exceeds", "TaoSetMaximumIterations", tao->max_it, &tao->max_it, &flg));
499: if (flg) PetscCall(TaoSetMaximumIterations(tao, tao->max_it));
501: PetscCall(PetscOptionsInt("-tao_max_funcs", "Stop if number of function evaluations exceeds", "TaoSetMaximumFunctionEvaluations", tao->max_funcs, &tao->max_funcs, &flg));
502: if (flg) PetscCall(TaoSetMaximumFunctionEvaluations(tao, tao->max_funcs));
504: PetscCall(PetscOptionsReal("-tao_fmin", "Stop if function less than", "TaoSetFunctionLowerBound", tao->fmin, &tao->fmin, NULL));
505: PetscCall(PetscOptionsBoundedReal("-tao_steptol", "Stop if step size or trust region radius less than", "", tao->steptol, &tao->steptol, NULL, 0));
506: PetscCall(PetscOptionsReal("-tao_trust0", "Initial trust region radius", "TaoSetInitialTrustRegionRadius", tao->trust0, &tao->trust0, &flg));
507: if (flg) PetscCall(TaoSetInitialTrustRegionRadius(tao, tao->trust0));
509: PetscCall(PetscOptionsDeprecated("-tao_solution_monitor", "-tao_monitor_solution", "3.21", NULL));
510: PetscCall(PetscOptionsDeprecated("-tao_gradient_monitor", "-tao_monitor_gradient", "3.21", NULL));
511: PetscCall(PetscOptionsDeprecated("-tao_stepdirection_monitor", "-tao_monitor_step", "3.21", NULL));
512: PetscCall(PetscOptionsDeprecated("-tao_residual_monitor", "-tao_monitor_residual", "3.21", NULL));
513: PetscCall(PetscOptionsDeprecated("-tao_smonitor", "-tao_monitor_short", "3.21", NULL));
514: PetscCall(PetscOptionsDeprecated("-tao_cmonitor", "-tao_monitor_constraint_norm", "3.21", NULL));
515: PetscCall(PetscOptionsDeprecated("-tao_gmonitor", "-tao_monitor_globalization", "3.21", NULL));
516: PetscCall(PetscOptionsDeprecated("-tao_draw_solution", "-tao_monitor_solution_draw", "3.21", NULL));
517: PetscCall(PetscOptionsDeprecated("-tao_draw_gradient", "-tao_monitor_gradient_draw", "3.21", NULL));
518: PetscCall(PetscOptionsDeprecated("-tao_draw_step", "-tao_monitor_step_draw", "3.21", NULL));
520: PetscCall(PetscOptionsBool("-tao_converged_reason", "Print reason for Tao converged", "TaoSolve", tao->printreason, &tao->printreason, NULL));
522: PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_solution", "View solution vector after each iteration", "TaoMonitorSolution", TaoMonitorSolution));
523: PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_gradient", "View gradient vector for each iteration", "TaoMonitorGradient", TaoMonitorGradient));
525: PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_step", "View step vector after each iteration", "TaoMonitorStep", TaoMonitorStep));
526: PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_residual", "View least-squares residual vector after each iteration", "TaoMonitorResidual", TaoMonitorResidual));
527: PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor", "Use the default convergence monitor", "TaoMonitorDefault", TaoMonitorDefault));
528: PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_globalization", "Use the convergence monitor with extra globalization info", "TaoMonitorGlobalization", TaoMonitorGlobalization));
529: PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_short", "Use the short convergence monitor", "TaoMonitorDefaultShort", TaoMonitorDefaultShort));
530: PetscCall(TaoMonitorSetFromOptions(tao, "-tao_monitor_constraint_norm", "Use the default convergence monitor with constraint norm", "TaoMonitorConstraintNorm", TaoMonitorConstraintNorm));
532: flg = PETSC_FALSE;
533: PetscCall(PetscOptionsDeprecated("-tao_cancelmonitors", "-tao_monitor_cancel", "3.21", NULL));
534: PetscCall(PetscOptionsBool("-tao_monitor_cancel", "cancel all monitors and call any registered destroy routines", "TaoMonitorCancel", flg, &flg, NULL));
535: if (flg) PetscCall(TaoMonitorCancel(tao));
537: flg = PETSC_FALSE;
538: PetscCall(PetscOptionsBool("-tao_monitor_solution_draw", "Plot solution vector at each iteration", "TaoMonitorSet", flg, &flg, NULL));
539: if (flg) {
540: TaoMonitorDrawCtx drawctx;
541: PetscInt howoften = 1;
542: PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &drawctx));
543: PetscCall(TaoMonitorSet(tao, TaoMonitorSolutionDraw, drawctx, (PetscCtxDestroyFn *)TaoMonitorDrawCtxDestroy));
544: }
546: flg = PETSC_FALSE;
547: PetscCall(PetscOptionsBool("-tao_monitor_step_draw", "Plots step at each iteration", "TaoMonitorSet", flg, &flg, NULL));
548: if (flg) PetscCall(TaoMonitorSet(tao, TaoMonitorStepDraw, NULL, NULL));
550: flg = PETSC_FALSE;
551: PetscCall(PetscOptionsBool("-tao_monitor_gradient_draw", "plots gradient at each iteration", "TaoMonitorSet", flg, &flg, NULL));
552: if (flg) {
553: TaoMonitorDrawCtx drawctx;
554: PetscInt howoften = 1;
555: PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &drawctx));
556: PetscCall(TaoMonitorSet(tao, TaoMonitorGradientDraw, drawctx, (PetscCtxDestroyFn *)TaoMonitorDrawCtxDestroy));
557: }
559: flg = PETSC_FALSE;
560: PetscCall(PetscOptionsBool("-tao_fd_gradient", "compute gradient using finite differences", "TaoDefaultComputeGradient", flg, &flg, NULL));
561: if (flg) PetscCall(TaoTermComputeGradientSetUseFD(tao->objective_term.term, PETSC_TRUE));
562: flg = PETSC_FALSE;
563: PetscCall(PetscOptionsBool("-tao_fd_hessian", "compute Hessian using finite differences", "TaoDefaultComputeHessian", flg, &flg, NULL));
564: if (flg) {
565: Mat H;
567: PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &H));
568: PetscCall(MatSetType(H, MATAIJ));
569: PetscCall(MatSetOption(H, MAT_SYMMETRIC, PETSC_TRUE));
570: PetscCall(MatSetOption(H, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
571: PetscCall(TaoSetHessian(tao, H, H, TaoDefaultComputeHessian, NULL));
572: PetscCall(TaoTermComputeHessianSetUseFD(tao->objective_term.term, PETSC_TRUE));
573: PetscCall(MatDestroy(&H));
574: }
575: flg = PETSC_FALSE;
576: PetscCall(PetscOptionsBool("-tao_mf_hessian", "compute matrix-free Hessian using finite differences", "TaoDefaultComputeHessianMFFD", flg, &flg, NULL));
577: if (flg) {
578: PetscBool is_callback;
579: Mat H;
581: // Check that tao has only one TaoTerm with type TAOTERMCALLBACK
582: PetscCall(PetscObjectTypeCompare((PetscObject)tao->objective_term.term, TAOTERMCALLBACKS, &is_callback));
583: if (is_callback) {
584: // Create Hessian via TaoTermCreateHessianMFFD
585: PetscCall(TaoTermCreateHessianMFFD(tao->objective_term.term, &H));
586: PetscCall(TaoSetHessian(tao, H, H, TaoDefaultComputeHessianMFFD, NULL));
587: PetscCall(MatDestroy(&H));
588: } else {
589: PetscCall(PetscInfo(tao, "-tao_mf_hessian only works when Tao has a single TAOTERMCALLBACK term. Ignoring.\n"));
590: }
591: }
592: PetscCall(PetscOptionsBool("-tao_recycle_history", "enable recycling/re-using information from the previous TaoSolve() call for some algorithms", "TaoSetRecycleHistory", flg, &flg, &found));
593: if (found) PetscCall(TaoSetRecycleHistory(tao, flg));
594: PetscCall(PetscOptionsEnum("-tao_subset_type", "subset type", "", TaoSubSetTypes, (PetscEnum)tao->subset_type, (PetscEnum *)&tao->subset_type, NULL));
596: if (tao->ksp) {
597: PetscCall(PetscOptionsBool("-tao_ksp_ew", "Use Eisentat-Walker linear system convergence test", "TaoKSPSetUseEW", tao->ksp_ewconv, &tao->ksp_ewconv, NULL));
598: PetscCall(TaoKSPSetUseEW(tao, tao->ksp_ewconv));
599: }
601: PetscCall(TaoTermSetFromOptions(tao->callbacks));
603: {
604: char *term_prefixes[16];
605: PetscInt n_terms = PETSC_STATIC_ARRAY_LENGTH(term_prefixes);
607: PetscCall(PetscOptionsStringArray("-tao_add_terms", "a list of prefixes for terms to add to the Tao objective function", "TaoAddTerm", term_prefixes, &n_terms, NULL));
608: for (PetscInt i = 0; i < n_terms; i++) {
609: TaoTerm term;
610: const char *prefix;
612: PetscCall(TaoTermDuplicate(tao->objective_term.term, TAOTERM_DUPLICATE_SIZEONLY, &term));
613: PetscCall(TaoGetOptionsPrefix(tao, &prefix));
614: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)term, prefix));
615: PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)term, term_prefixes[i]));
616: PetscCall(TaoTermSetFromOptions(term));
617: PetscCall(TaoAddTerm(tao, term_prefixes[i], 1.0, term, NULL, NULL));
618: PetscCall(TaoTermDestroy(&term));
619: PetscCall(PetscFree(term_prefixes[i]));
620: }
621: }
623: if (tao->objective_term.term != tao->callbacks) PetscCall(TaoTermSetFromOptions(tao->objective_term.term));
625: PetscTryTypeMethod(tao, setfromoptions, PetscOptionsObject);
627: /* process any options handlers added with PetscObjectAddOptionsHandler() */
628: PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)tao, PetscOptionsObject));
629: PetscOptionsEnd();
631: if (tao->linesearch) PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
632: PetscFunctionReturn(PETSC_SUCCESS);
633: }
635: /*@
636: TaoViewFromOptions - View a `Tao` object based on values in the options database
638: Collective
640: Input Parameters:
641: + A - the `Tao` context
642: . obj - Optional object that provides the prefix for the options database
643: - name - command line option
645: Options Database Key:
646: . -name [viewertype][:...] - option name and values. See `PetscObjectViewFromOptions()` for the possible arguments
648: Level: intermediate
650: .seealso: [](ch_tao), `Tao`, `TaoView`, `PetscObjectViewFromOptions()`, `TaoCreate()`
651: @*/
652: PetscErrorCode TaoViewFromOptions(Tao A, PetscObject obj, const char name[])
653: {
654: PetscFunctionBegin;
656: PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
657: PetscFunctionReturn(PETSC_SUCCESS);
658: }
660: /*@
661: TaoView - Prints information about the `Tao` object
663: Collective
665: Input Parameters:
666: + tao - the `Tao` context
667: - viewer - visualization context
669: Options Database Key:
670: . -tao_view - Calls `TaoView()` at the end of `TaoSolve()`
672: Level: beginner
674: Notes:
675: The available visualization contexts include
676: + `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
677: - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
678: output where only the first processor opens
679: the file. All other processors send their
680: data to the first processor to print.
682: To view all the `TaoTerm` inside of `Tao`, use `PETSC_VIEWER_ASCII_INFO_DETAIL`,
683: or pass `-tao_view ::ascii_info_detail` flag
685: .seealso: [](ch_tao), `Tao`, `PetscViewerASCIIOpen()`
686: @*/
687: PetscErrorCode TaoView(Tao tao, PetscViewer viewer)
688: {
689: PetscBool isascii, isstring;
690: TaoType type;
692: PetscFunctionBegin;
694: if (!viewer) PetscCall(PetscViewerASCIIGetStdout(((PetscObject)tao)->comm, &viewer));
696: PetscCheckSameComm(tao, 1, viewer, 2);
698: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
699: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
700: if (isascii) {
701: PetscViewerFormat format;
703: PetscCall(PetscViewerGetFormat(viewer, &format));
704: PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)tao, viewer));
706: PetscCall(PetscViewerASCIIPushTab(viewer));
707: PetscTryTypeMethod(tao, view, viewer);
708: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
709: PetscCall(PetscViewerASCIIPrintf(viewer, "Objective function:\n"));
710: PetscCall(PetscViewerASCIIPushTab(viewer));
711: PetscCall(PetscViewerASCIIPrintf(viewer, "Scale (tao_objective_scale): %g\n", (double)tao->objective_term.scale));
712: PetscCall(PetscViewerASCIIPrintf(viewer, "Function:\n"));
713: PetscCall(PetscViewerASCIIPushTab(viewer));
714: PetscCall(TaoTermView(tao->objective_term.term, viewer));
715: PetscCall(PetscViewerASCIIPopTab(viewer));
716: if (tao->objective_term.map) {
717: PetscCall(PetscViewerASCIIPrintf(viewer, "Map:\n"));
718: PetscCall(PetscViewerASCIIPushTab(viewer));
719: PetscCall(MatView(tao->objective_term.map, viewer));
720: PetscCall(PetscViewerASCIIPopTab(viewer));
721: } else PetscCall(PetscViewerASCIIPrintf(viewer, "Map: unmapped\n"));
722: PetscCall(PetscViewerASCIIPopTab(viewer));
723: } else if (tao->num_terms > 0 || tao->term_set) {
724: if (tao->objective_term.scale == 1.0 && tao->objective_term.map == NULL) {
725: PetscCall(PetscViewerASCIIPrintf(viewer, "Objective function:\n"));
726: PetscCall(PetscViewerASCIIPushTab(viewer));
727: PetscCall(TaoTermView(tao->objective_term.term, viewer));
728: PetscCall(PetscViewerASCIIPopTab(viewer));
729: } else {
730: PetscCall(PetscViewerASCIIPrintf(viewer, "Objective function:\n"));
731: PetscCall(PetscViewerASCIIPushTab(viewer));
732: if (tao->objective_term.scale != 1.0) PetscCall(PetscViewerASCIIPrintf(viewer, "Scale: %g\n", (double)tao->objective_term.scale));
733: PetscCall(PetscViewerASCIIPrintf(viewer, "Function:\n"));
734: PetscCall(PetscViewerASCIIPushTab(viewer));
735: PetscCall(TaoTermView(tao->objective_term.term, viewer));
736: PetscCall(PetscViewerASCIIPopTab(viewer));
737: if (tao->objective_term.map) {
738: PetscCall(PetscViewerASCIIPrintf(viewer, "Map:\n"));
739: PetscCall(PetscViewerASCIIPushTab(viewer));
740: PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO));
741: PetscCall(MatView(tao->objective_term.map, viewer));
742: PetscCall(PetscViewerPopFormat(viewer));
743: PetscCall(PetscViewerASCIIPopTab(viewer));
744: }
745: PetscCall(PetscViewerASCIIPopTab(viewer));
746: }
747: }
748: if (tao->linesearch) PetscCall(TaoLineSearchView(tao->linesearch, viewer));
749: if (tao->ksp) {
750: PetscCall(KSPView(tao->ksp, viewer));
751: PetscCall(PetscViewerASCIIPrintf(viewer, "total KSP iterations: %" PetscInt_FMT "\n", tao->ksp_tot_its));
752: }
754: if (tao->XL || tao->XU) PetscCall(PetscViewerASCIIPrintf(viewer, "Active Set subset type: %s\n", TaoSubSetTypes[tao->subset_type]));
756: PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: gatol=%g,", (double)tao->gatol));
757: PetscCall(PetscViewerASCIIPrintf(viewer, " grtol=%g,", (double)tao->grtol));
758: PetscCall(PetscViewerASCIIPrintf(viewer, " steptol=%g,", (double)tao->steptol));
759: PetscCall(PetscViewerASCIIPrintf(viewer, " gttol=%g\n", (double)tao->gttol));
760: PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Function/Gradient:=%g\n", (double)tao->residual));
762: if (tao->constrained) {
763: PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances:"));
764: PetscCall(PetscViewerASCIIPrintf(viewer, " catol=%g,", (double)tao->catol));
765: PetscCall(PetscViewerASCIIPrintf(viewer, " crtol=%g\n", (double)tao->crtol));
766: PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Constraints:=%g\n", (double)tao->cnorm));
767: }
769: if (tao->trust < tao->steptol) {
770: PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: steptol=%g\n", (double)tao->steptol));
771: PetscCall(PetscViewerASCIIPrintf(viewer, "Final trust region radius:=%g\n", (double)tao->trust));
772: }
774: if (tao->fmin > -1.e25) PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: function minimum=%g\n", (double)tao->fmin));
775: PetscCall(PetscViewerASCIIPrintf(viewer, "Objective value=%g\n", (double)tao->fc));
777: PetscCall(PetscViewerASCIIPrintf(viewer, "total number of iterations=%" PetscInt_FMT ", ", tao->niter));
778: PetscCall(PetscViewerASCIIPrintf(viewer, " (max: %" PetscInt_FMT ")\n", tao->max_it));
780: if (tao->objective_term.term->nobj > 0) {
781: PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function evaluations=%" PetscInt_FMT ",", tao->objective_term.term->nobj));
782: if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, " (max: unlimited)\n"));
783: else PetscCall(PetscViewerASCIIPrintf(viewer, " (max: %" PetscInt_FMT ")\n", tao->max_funcs));
784: }
785: if (tao->objective_term.term->ngrad > 0) {
786: PetscCall(PetscViewerASCIIPrintf(viewer, "total number of gradient evaluations=%" PetscInt_FMT ",", tao->objective_term.term->ngrad));
787: if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, " (max: unlimited)\n"));
788: else PetscCall(PetscViewerASCIIPrintf(viewer, " (max: %" PetscInt_FMT ")\n", tao->max_funcs));
789: }
790: if (tao->objective_term.term->nobjgrad > 0) {
791: PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function/gradient evaluations=%" PetscInt_FMT ",", tao->objective_term.term->nobjgrad));
792: if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, " (max: unlimited)\n"));
793: else PetscCall(PetscViewerASCIIPrintf(viewer, " (max: %" PetscInt_FMT ")\n", tao->max_funcs));
794: }
795: if (tao->nres > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of residual evaluations=%" PetscInt_FMT "\n", tao->nres));
796: if (tao->objective_term.term->nhess > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Hessian evaluations=%" PetscInt_FMT "\n", tao->objective_term.term->nhess));
797: if (tao->nconstraints > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of constraint function evaluations=%" PetscInt_FMT "\n", tao->nconstraints));
798: if (tao->njac > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Jacobian evaluations=%" PetscInt_FMT "\n", tao->njac));
800: if (tao->reason > 0) {
801: PetscCall(PetscViewerASCIIPrintf(viewer, "Solution converged: "));
802: switch (tao->reason) {
803: case TAO_CONVERGED_GATOL:
804: PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)|| <= gatol\n"));
805: break;
806: case TAO_CONVERGED_GRTOL:
807: PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)||/|f(X)| <= grtol\n"));
808: break;
809: case TAO_CONVERGED_GTTOL:
810: PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)||/||g(X0)|| <= gttol\n"));
811: break;
812: case TAO_CONVERGED_STEPTOL:
813: PetscCall(PetscViewerASCIIPrintf(viewer, " Steptol -- step size small\n"));
814: break;
815: case TAO_CONVERGED_MINF:
816: PetscCall(PetscViewerASCIIPrintf(viewer, " Minf -- f < fmin\n"));
817: break;
818: case TAO_CONVERGED_USER:
819: PetscCall(PetscViewerASCIIPrintf(viewer, " User Terminated\n"));
820: break;
821: default:
822: PetscCall(PetscViewerASCIIPrintf(viewer, " %d\n", tao->reason));
823: break;
824: }
825: } else if (tao->reason == TAO_CONTINUE_ITERATING) {
826: PetscCall(PetscViewerASCIIPrintf(viewer, "Solver never run\n"));
827: } else {
828: PetscCall(PetscViewerASCIIPrintf(viewer, "Solver failed: "));
829: switch (tao->reason) {
830: case TAO_DIVERGED_MAXITS:
831: PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum Iterations\n"));
832: break;
833: case TAO_DIVERGED_NAN:
834: PetscCall(PetscViewerASCIIPrintf(viewer, " NAN or infinity encountered\n"));
835: break;
836: case TAO_DIVERGED_MAXFCN:
837: PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum Function Evaluations\n"));
838: break;
839: case TAO_DIVERGED_LS_FAILURE:
840: PetscCall(PetscViewerASCIIPrintf(viewer, " Line Search Failure\n"));
841: break;
842: case TAO_DIVERGED_TR_REDUCTION:
843: PetscCall(PetscViewerASCIIPrintf(viewer, " Trust Region too small\n"));
844: break;
845: case TAO_DIVERGED_USER:
846: PetscCall(PetscViewerASCIIPrintf(viewer, " User Terminated\n"));
847: break;
848: default:
849: PetscCall(PetscViewerASCIIPrintf(viewer, " %d\n", tao->reason));
850: break;
851: }
852: }
853: PetscCall(PetscViewerASCIIPopTab(viewer));
854: } else if (isstring) {
855: PetscCall(TaoGetType(tao, &type));
856: PetscCall(PetscViewerStringSPrintf(viewer, " %-3.3s", type));
857: }
858: PetscFunctionReturn(PETSC_SUCCESS);
859: }
861: /*@
862: TaoSetRecycleHistory - Sets the boolean flag to enable/disable re-using
863: iterate information from the previous `TaoSolve()`. This feature is disabled by
864: default.
866: Logically Collective
868: Input Parameters:
869: + tao - the `Tao` context
870: - recycle - boolean flag
872: Options Database Key:
873: . -tao_recycle_history (true|false) - reuse the history
875: Level: intermediate
877: Notes:
878: For conjugate gradient methods (`TAOBNCG`), this re-uses the latest search direction
879: from the previous `TaoSolve()` call when computing the first search direction in a
880: new solution. By default, CG methods set the first search direction to the
881: negative gradient.
883: For quasi-Newton family of methods (`TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`), this re-uses
884: the accumulated quasi-Newton Hessian approximation from the previous `TaoSolve()`
885: call. By default, QN family of methods reset the initial Hessian approximation to
886: the identity matrix.
888: For any other algorithm, this setting has no effect.
890: .seealso: [](ch_tao), `Tao`, `TaoGetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
891: @*/
892: PetscErrorCode TaoSetRecycleHistory(Tao tao, PetscBool recycle)
893: {
894: PetscFunctionBegin;
897: tao->recycle = recycle;
898: PetscFunctionReturn(PETSC_SUCCESS);
899: }
901: /*@
902: TaoGetRecycleHistory - Retrieve the boolean flag for re-using iterate information
903: from the previous `TaoSolve()`. This feature is disabled by default.
905: Logically Collective
907: Input Parameter:
908: . tao - the `Tao` context
910: Output Parameter:
911: . recycle - boolean flag
913: Level: intermediate
915: .seealso: [](ch_tao), `Tao`, `TaoSetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
916: @*/
917: PetscErrorCode TaoGetRecycleHistory(Tao tao, PetscBool *recycle)
918: {
919: PetscFunctionBegin;
921: PetscAssertPointer(recycle, 2);
922: *recycle = tao->recycle;
923: PetscFunctionReturn(PETSC_SUCCESS);
924: }
926: /*@
927: TaoSetTolerances - Sets parameters used in `TaoSolve()` convergence tests
929: Logically Collective
931: Input Parameters:
932: + tao - the `Tao` context
933: . gatol - stop if norm of gradient is less than this
934: . grtol - stop if relative norm of gradient is less than this
935: - gttol - stop if norm of gradient is reduced by this factor
937: Options Database Keys:
938: + -tao_gatol gatol - Sets gatol
939: . -tao_grtol grtol - Sets grtol
940: - -tao_gttol gttol - Sets gttol
942: Stopping Criteria\:
943: .vb
944: ||g(X)|| <= gatol
945: ||g(X)|| / |f(X)| <= grtol
946: ||g(X)|| / ||g(X0)|| <= gttol
947: .ve
949: Level: beginner
951: Notes:
952: Use `PETSC_CURRENT` to leave one or more tolerances unchanged.
954: Use `PETSC_DETERMINE` to set one or more tolerances to their values when the `tao`object's type was set
956: Fortran Note:
957: Use `PETSC_CURRENT_REAL` or `PETSC_DETERMINE_REAL`
959: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`
960: @*/
961: PetscErrorCode TaoSetTolerances(Tao tao, PetscReal gatol, PetscReal grtol, PetscReal gttol)
962: {
963: PetscFunctionBegin;
969: if (gatol == (PetscReal)PETSC_DETERMINE) {
970: tao->gatol = tao->default_gatol;
971: } else if (gatol != (PetscReal)PETSC_CURRENT) {
972: PetscCheck(gatol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative gatol not allowed");
973: tao->gatol = gatol;
974: }
976: if (grtol == (PetscReal)PETSC_DETERMINE) {
977: tao->grtol = tao->default_grtol;
978: } else if (grtol != (PetscReal)PETSC_CURRENT) {
979: PetscCheck(grtol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative grtol not allowed");
980: tao->grtol = grtol;
981: }
983: if (gttol == (PetscReal)PETSC_DETERMINE) {
984: tao->gttol = tao->default_gttol;
985: } else if (gttol != (PetscReal)PETSC_CURRENT) {
986: PetscCheck(gttol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative gttol not allowed");
987: tao->gttol = gttol;
988: }
989: PetscFunctionReturn(PETSC_SUCCESS);
990: }
992: /*@
993: TaoSetConstraintTolerances - Sets constraint tolerance parameters used in `TaoSolve()` convergence tests
995: Logically Collective
997: Input Parameters:
998: + tao - the `Tao` context
999: . catol - absolute constraint tolerance, constraint norm must be less than `catol` for used for `gatol` convergence criteria
1000: - crtol - relative constraint tolerance, constraint norm must be less than `crtol` for used for `gatol`, `gttol` convergence criteria
1002: Options Database Keys:
1003: + -tao_catol catol - Sets catol
1004: - -tao_crtol crtol - Sets crtol
1006: Level: intermediate
1008: Notes:
1009: Use `PETSC_CURRENT` to leave one or tolerance unchanged.
1011: Use `PETSC_DETERMINE` to set one or more tolerances to their values when the `tao` object's type was set
1013: Fortran Note:
1014: Use `PETSC_CURRENT_REAL` or `PETSC_DETERMINE_REAL`
1016: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`, `TaoGetConstraintTolerances()`, `TaoSetTolerances()`
1017: @*/
1018: PetscErrorCode TaoSetConstraintTolerances(Tao tao, PetscReal catol, PetscReal crtol)
1019: {
1020: PetscFunctionBegin;
1025: if (catol == (PetscReal)PETSC_DETERMINE) {
1026: tao->catol = tao->default_catol;
1027: } else if (catol != (PetscReal)PETSC_CURRENT) {
1028: PetscCheck(catol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative catol not allowed");
1029: tao->catol = catol;
1030: }
1032: if (crtol == (PetscReal)PETSC_DETERMINE) {
1033: tao->crtol = tao->default_crtol;
1034: } else if (crtol != (PetscReal)PETSC_CURRENT) {
1035: PetscCheck(crtol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative crtol not allowed");
1036: tao->crtol = crtol;
1037: }
1038: PetscFunctionReturn(PETSC_SUCCESS);
1039: }
1041: /*@
1042: TaoGetConstraintTolerances - Gets constraint tolerance parameters used in `TaoSolve()` convergence tests
1044: Not Collective
1046: Input Parameter:
1047: . tao - the `Tao` context
1049: Output Parameters:
1050: + catol - absolute constraint tolerance, constraint norm must be less than `catol` for used for `gatol` convergence criteria
1051: - crtol - relative constraint tolerance, constraint norm must be less than `crtol` for used for `gatol`, `gttol` convergence criteria
1053: Level: intermediate
1055: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`, `TaoSetTolerances()`, `TaoSetConstraintTolerances()`
1056: @*/
1057: PetscErrorCode TaoGetConstraintTolerances(Tao tao, PetscReal *catol, PetscReal *crtol)
1058: {
1059: PetscFunctionBegin;
1061: if (catol) *catol = tao->catol;
1062: if (crtol) *crtol = tao->crtol;
1063: PetscFunctionReturn(PETSC_SUCCESS);
1064: }
1066: /*@
1067: TaoSetFunctionLowerBound - Sets a bound on the solution objective value.
1068: When an approximate solution with an objective value below this number
1069: has been found, the solver will terminate.
1071: Logically Collective
1073: Input Parameters:
1074: + tao - the Tao solver context
1075: - fmin - the tolerance
1077: Options Database Key:
1078: . -tao_fmin fmin - sets the minimum function value
1080: Level: intermediate
1082: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetTolerances()`
1083: @*/
1084: PetscErrorCode TaoSetFunctionLowerBound(Tao tao, PetscReal fmin)
1085: {
1086: PetscFunctionBegin;
1089: tao->fmin = fmin;
1090: PetscFunctionReturn(PETSC_SUCCESS);
1091: }
1093: /*@
1094: TaoGetFunctionLowerBound - Gets the bound on the solution objective value.
1095: When an approximate solution with an objective value below this number
1096: has been found, the solver will terminate.
1098: Not Collective
1100: Input Parameter:
1101: . tao - the `Tao` solver context
1103: Output Parameter:
1104: . fmin - the minimum function value
1106: Level: intermediate
1108: .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetFunctionLowerBound()`
1109: @*/
1110: PetscErrorCode TaoGetFunctionLowerBound(Tao tao, PetscReal *fmin)
1111: {
1112: PetscFunctionBegin;
1114: PetscAssertPointer(fmin, 2);
1115: *fmin = tao->fmin;
1116: PetscFunctionReturn(PETSC_SUCCESS);
1117: }
1119: /*@
1120: TaoSetMaximumFunctionEvaluations - Sets a maximum number of function evaluations allowed for a `TaoSolve()`.
1122: Logically Collective
1124: Input Parameters:
1125: + tao - the `Tao` solver context
1126: - nfcn - the maximum number of function evaluations (>=0), use `PETSC_UNLIMITED` to have no bound
1128: Options Database Key:
1129: . -tao_max_funcs nfcn - sets the maximum number of function evaluations
1131: Level: intermediate
1133: Note:
1134: Use `PETSC_DETERMINE` to use the default maximum number of function evaluations that was set when the object type was set.
1136: Developer Note:
1137: Deprecated support for an unlimited number of function evaluations by passing a negative value.
1139: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoSetMaximumIterations()`
1140: @*/
1141: PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao tao, PetscInt nfcn)
1142: {
1143: PetscFunctionBegin;
1146: if (nfcn == PETSC_DETERMINE) {
1147: tao->max_funcs = tao->default_max_funcs;
1148: } else if (nfcn == PETSC_UNLIMITED || nfcn < 0) {
1149: tao->max_funcs = PETSC_UNLIMITED;
1150: } else {
1151: PetscCheck(nfcn >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of function evaluations must be positive");
1152: tao->max_funcs = nfcn;
1153: }
1154: PetscFunctionReturn(PETSC_SUCCESS);
1155: }
1157: /*@
1158: TaoGetMaximumFunctionEvaluations - Gets a maximum number of function evaluations allowed for a `TaoSolve()`
1160: Logically Collective
1162: Input Parameter:
1163: . tao - the `Tao` solver context
1165: Output Parameter:
1166: . nfcn - the maximum number of function evaluations
1168: Level: intermediate
1170: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1171: @*/
1172: PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao tao, PetscInt *nfcn)
1173: {
1174: PetscFunctionBegin;
1176: PetscAssertPointer(nfcn, 2);
1177: *nfcn = tao->max_funcs;
1178: PetscFunctionReturn(PETSC_SUCCESS);
1179: }
1181: /*@
1182: TaoGetCurrentFunctionEvaluations - Get current number of function evaluations used by a `Tao` object
1184: Not Collective
1186: Input Parameter:
1187: . tao - the `Tao` solver context
1189: Output Parameter:
1190: . nfuncs - the current number of function evaluations (maximum between gradient and function evaluations)
1192: Level: intermediate
1194: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1195: @*/
1196: PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao tao, PetscInt *nfuncs)
1197: {
1198: PetscFunctionBegin;
1200: PetscAssertPointer(nfuncs, 2);
1201: *nfuncs = PetscMax(tao->objective_term.term->nobj, tao->objective_term.term->nobjgrad);
1202: PetscFunctionReturn(PETSC_SUCCESS);
1203: }
1205: /*@
1206: TaoSetMaximumIterations - Sets a maximum number of iterates to be used in `TaoSolve()`
1208: Logically Collective
1210: Input Parameters:
1211: + tao - the `Tao` solver context
1212: - maxits - the maximum number of iterates (>=0), use `PETSC_UNLIMITED` to have no bound
1214: Options Database Key:
1215: . -tao_max_it its - sets the maximum number of iterations
1217: Level: intermediate
1219: Note:
1220: Use `PETSC_DETERMINE` to use the default maximum number of iterations that was set when the object's type was set.
1222: Developer Note:
1223: Also accepts the deprecated negative values to indicate no limit
1225: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoSetMaximumFunctionEvaluations()`
1226: @*/
1227: PetscErrorCode TaoSetMaximumIterations(Tao tao, PetscInt maxits)
1228: {
1229: PetscFunctionBegin;
1232: if (maxits == PETSC_DETERMINE) {
1233: tao->max_it = tao->default_max_it;
1234: } else if (maxits == PETSC_UNLIMITED) {
1235: tao->max_it = PETSC_INT_MAX;
1236: } else {
1237: PetscCheck(maxits > 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations must be positive");
1238: tao->max_it = maxits;
1239: }
1240: PetscFunctionReturn(PETSC_SUCCESS);
1241: }
1243: /*@
1244: TaoGetMaximumIterations - Gets a maximum number of iterates that will be used
1246: Not Collective
1248: Input Parameter:
1249: . tao - the `Tao` solver context
1251: Output Parameter:
1252: . maxits - the maximum number of iterates
1254: Level: intermediate
1256: .seealso: [](ch_tao), `Tao`, `TaoSetMaximumIterations()`, `TaoGetMaximumFunctionEvaluations()`
1257: @*/
1258: PetscErrorCode TaoGetMaximumIterations(Tao tao, PetscInt *maxits)
1259: {
1260: PetscFunctionBegin;
1262: PetscAssertPointer(maxits, 2);
1263: *maxits = tao->max_it;
1264: PetscFunctionReturn(PETSC_SUCCESS);
1265: }
1267: /*@
1268: TaoSetInitialTrustRegionRadius - Sets the initial trust region radius.
1270: Logically Collective
1272: Input Parameters:
1273: + tao - a `Tao` optimization solver
1274: - radius - the trust region radius
1276: Options Database Key:
1277: . -tao_trust0 radius - sets initial trust region radius
1279: Level: intermediate
1281: Note:
1282: Use `PETSC_DETERMINE` to use the default radius that was set when the object's type was set.
1284: .seealso: [](ch_tao), `Tao`, `TaoGetTrustRegionRadius()`, `TaoSetTrustRegionTolerance()`, `TAONTR`
1285: @*/
1286: PetscErrorCode TaoSetInitialTrustRegionRadius(Tao tao, PetscReal radius)
1287: {
1288: PetscFunctionBegin;
1291: if (radius == PETSC_DETERMINE) {
1292: tao->trust0 = tao->default_trust0;
1293: } else {
1294: PetscCheck(radius > 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Radius must be positive");
1295: tao->trust0 = radius;
1296: }
1297: PetscFunctionReturn(PETSC_SUCCESS);
1298: }
1300: /*@
1301: TaoGetInitialTrustRegionRadius - Gets the initial trust region radius.
1303: Not Collective
1305: Input Parameter:
1306: . tao - a `Tao` optimization solver
1308: Output Parameter:
1309: . radius - the trust region radius
1311: Level: intermediate
1313: .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetCurrentTrustRegionRadius()`, `TAONTR`
1314: @*/
1315: PetscErrorCode TaoGetInitialTrustRegionRadius(Tao tao, PetscReal *radius)
1316: {
1317: PetscFunctionBegin;
1319: PetscAssertPointer(radius, 2);
1320: *radius = tao->trust0;
1321: PetscFunctionReturn(PETSC_SUCCESS);
1322: }
1324: /*@
1325: TaoGetCurrentTrustRegionRadius - Gets the current trust region radius.
1327: Not Collective
1329: Input Parameter:
1330: . tao - a `Tao` optimization solver
1332: Output Parameter:
1333: . radius - the trust region radius
1335: Level: intermediate
1337: .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetInitialTrustRegionRadius()`, `TAONTR`
1338: @*/
1339: PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao tao, PetscReal *radius)
1340: {
1341: PetscFunctionBegin;
1343: PetscAssertPointer(radius, 2);
1344: *radius = tao->trust;
1345: PetscFunctionReturn(PETSC_SUCCESS);
1346: }
1348: /*@
1349: TaoGetTolerances - gets the current values of some tolerances used for the convergence testing of `TaoSolve()`
1351: Not Collective
1353: Input Parameter:
1354: . tao - the `Tao` context
1356: Output Parameters:
1357: + gatol - stop if norm of gradient is less than this
1358: . grtol - stop if relative norm of gradient is less than this
1359: - gttol - stop if norm of gradient is reduced by a this factor
1361: Level: intermediate
1363: Note:
1364: `NULL` can be used as an argument if not all tolerances values are needed
1366: .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`
1367: @*/
1368: PetscErrorCode TaoGetTolerances(Tao tao, PetscReal *gatol, PetscReal *grtol, PetscReal *gttol)
1369: {
1370: PetscFunctionBegin;
1372: if (gatol) *gatol = tao->gatol;
1373: if (grtol) *grtol = tao->grtol;
1374: if (gttol) *gttol = tao->gttol;
1375: PetscFunctionReturn(PETSC_SUCCESS);
1376: }
1378: /*@
1379: TaoGetKSP - Gets the linear solver used by the optimization solver.
1381: Not Collective
1383: Input Parameter:
1384: . tao - the `Tao` solver
1386: Output Parameter:
1387: . ksp - the `KSP` linear solver used in the optimization solver
1389: Level: intermediate
1391: .seealso: [](ch_tao), `Tao`, `KSP`
1392: @*/
1393: PetscErrorCode TaoGetKSP(Tao tao, KSP *ksp)
1394: {
1395: PetscFunctionBegin;
1397: PetscAssertPointer(ksp, 2);
1398: *ksp = tao->ksp;
1399: PetscFunctionReturn(PETSC_SUCCESS);
1400: }
1402: /*@
1403: TaoGetLinearSolveIterations - Gets the total number of linear iterations
1404: used by the `Tao` solver
1406: Not Collective
1408: Input Parameter:
1409: . tao - the `Tao` context
1411: Output Parameter:
1412: . lits - number of linear iterations
1414: Level: intermediate
1416: Note:
1417: This counter is reset to zero for each successive call to `TaoSolve()`
1419: .seealso: [](ch_tao), `Tao`, `TaoGetKSP()`
1420: @*/
1421: PetscErrorCode TaoGetLinearSolveIterations(Tao tao, PetscInt *lits)
1422: {
1423: PetscFunctionBegin;
1425: PetscAssertPointer(lits, 2);
1426: *lits = tao->ksp_tot_its;
1427: PetscFunctionReturn(PETSC_SUCCESS);
1428: }
1430: /*@
1431: TaoGetLineSearch - Gets the line search used by the optimization solver.
1433: Not Collective
1435: Input Parameter:
1436: . tao - the `Tao` solver
1438: Output Parameter:
1439: . ls - the line search used in the optimization solver
1441: Level: intermediate
1443: .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`
1444: @*/
1445: PetscErrorCode TaoGetLineSearch(Tao tao, TaoLineSearch *ls)
1446: {
1447: PetscFunctionBegin;
1449: PetscAssertPointer(ls, 2);
1450: *ls = tao->linesearch;
1451: PetscFunctionReturn(PETSC_SUCCESS);
1452: }
1454: /*@
1455: TaoAddLineSearchCounts - Adds the number of function evaluations spent
1456: in the line search to the running total.
1458: Input Parameters:
1459: . tao - the `Tao` solver
1461: Level: developer
1463: .seealso: [](ch_tao), `Tao`, `TaoGetLineSearch()`, `TaoLineSearchApply()`
1464: @*/
1465: PetscErrorCode TaoAddLineSearchCounts(Tao tao)
1466: {
1467: PetscBool flg;
1468: PetscInt nfeval, ngeval, nfgeval;
1470: PetscFunctionBegin;
1472: if (tao->linesearch) {
1473: PetscCall(TaoLineSearchIsUsingTaoRoutines(tao->linesearch, &flg));
1474: if (!flg) {
1475: PetscCall(TaoLineSearchGetNumberFunctionEvaluations(tao->linesearch, &nfeval, &ngeval, &nfgeval));
1476: tao->objective_term.term->nobj += nfeval;
1477: tao->objective_term.term->ngrad += ngeval;
1478: tao->objective_term.term->nobjgrad += nfgeval;
1479: }
1480: }
1481: PetscFunctionReturn(PETSC_SUCCESS);
1482: }
1484: /*@
1485: TaoGetSolution - Returns the vector with the current solution from the `Tao` object
1487: Not Collective
1489: Input Parameter:
1490: . tao - the `Tao` context
1492: Output Parameter:
1493: . X - the current solution
1495: Level: intermediate
1497: Note:
1498: The returned vector will be the same object that was passed into `TaoSetSolution()`
1500: .seealso: [](ch_tao), `Tao`, `TaoSetSolution()`, `TaoSolve()`
1501: @*/
1502: PetscErrorCode TaoGetSolution(Tao tao, Vec *X)
1503: {
1504: PetscFunctionBegin;
1506: PetscAssertPointer(X, 2);
1507: *X = tao->solution;
1508: PetscFunctionReturn(PETSC_SUCCESS);
1509: }
1511: /*@
1512: TaoResetStatistics - Initialize the statistics collected by the `Tao` object.
1513: These statistics include the iteration number, residual norms, and convergence status.
1514: This routine gets called before solving each optimization problem.
1516: Collective
1518: Input Parameter:
1519: . tao - the `Tao` context
1521: Level: developer
1523: Note:
1524: This function does not reset the statistics of internal `TaoTerm`
1526: .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
1527: @*/
1528: PetscErrorCode TaoResetStatistics(Tao tao)
1529: {
1530: PetscFunctionBegin;
1532: tao->niter = 0;
1533: tao->nres = 0;
1534: tao->njac = 0;
1535: tao->nconstraints = 0;
1536: tao->ksp_its = 0;
1537: tao->ksp_tot_its = 0;
1538: tao->reason = TAO_CONTINUE_ITERATING;
1539: tao->residual = 0.0;
1540: tao->cnorm = 0.0;
1541: tao->step = 0.0;
1542: tao->lsflag = PETSC_FALSE;
1543: if (tao->hist_reset) tao->hist_len = 0;
1544: PetscFunctionReturn(PETSC_SUCCESS);
1545: }
1547: /*@C
1548: TaoSetUpdate - Sets the general-purpose update function called
1549: at the beginning of every iteration of the optimization algorithm. Called after the new solution and the gradient
1550: is determined, but before the Hessian is computed (if applicable).
1552: Logically Collective
1554: Input Parameters:
1555: + tao - The `Tao` solver
1556: . func - The function
1557: - ctx - The update function context
1559: Calling sequence of `func`:
1560: + tao - The optimizer context
1561: . it - The current iteration index
1562: - ctx - The update context
1564: Level: advanced
1566: Notes:
1567: Users can modify the gradient direction or any other vector associated to the specific solver used.
1568: The objective function value is always recomputed after a call to the update hook.
1570: .seealso: [](ch_tao), `Tao`, `TaoSolve()`
1571: @*/
1572: PetscErrorCode TaoSetUpdate(Tao tao, PetscErrorCode (*func)(Tao tao, PetscInt it, PetscCtx ctx), PetscCtx ctx)
1573: {
1574: PetscFunctionBegin;
1576: tao->ops->update = func;
1577: tao->user_update = ctx;
1578: PetscFunctionReturn(PETSC_SUCCESS);
1579: }
1581: /*@C
1582: TaoSetConvergenceTest - Sets the function that is to be used to test
1583: for convergence of the iterative minimization solution. The new convergence
1584: testing routine will replace Tao's default convergence test.
1586: Logically Collective
1588: Input Parameters:
1589: + tao - the `Tao` object
1590: . conv - the routine to test for convergence
1591: - ctx - [optional] context for private data for the convergence routine (may be `NULL`)
1593: Calling sequence of `conv`:
1594: + tao - the `Tao` object
1595: - ctx - [optional] convergence context
1597: Level: advanced
1599: Note:
1600: The new convergence testing routine should call `TaoSetConvergedReason()`.
1602: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetConvergedReason()`, `TaoGetSolutionStatus()`, `TaoGetTolerances()`, `TaoMonitorSet()`
1603: @*/
1604: PetscErrorCode TaoSetConvergenceTest(Tao tao, PetscErrorCode (*conv)(Tao tao, PetscCtx ctx), PetscCtx ctx)
1605: {
1606: PetscFunctionBegin;
1608: tao->ops->convergencetest = conv;
1609: tao->cnvP = ctx;
1610: PetscFunctionReturn(PETSC_SUCCESS);
1611: }
1613: /*@C
1614: TaoMonitorSet - Sets an additional function that is to be used at every
1615: iteration of the solver to display the iteration's
1616: progress.
1618: Logically Collective
1620: Input Parameters:
1621: + tao - the `Tao` solver context
1622: . func - monitoring routine
1623: . ctx - [optional] user-defined context for private data for the monitor routine (may be `NULL`)
1624: - dest - [optional] function to destroy the context when the `Tao` is destroyed, see `PetscCtxDestroyFn` for the calling sequence
1626: Calling sequence of `func`:
1627: + tao - the `Tao` solver context
1628: - ctx - [optional] monitoring context
1630: Level: intermediate
1632: Notes:
1633: See `TaoSetFromOptions()` for a monitoring options.
1635: Several different monitoring routines may be set by calling
1636: `TaoMonitorSet()` multiple times; all will be called in the
1637: order in which they were set.
1639: Fortran Notes:
1640: Only one monitor function may be set
1642: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoMonitorDefault()`, `TaoMonitorCancel()`, `TaoView()`, `PetscCtxDestroyFn`
1643: @*/
1644: PetscErrorCode TaoMonitorSet(Tao tao, PetscErrorCode (*func)(Tao tao, PetscCtx ctx), PetscCtx ctx, PetscCtxDestroyFn *dest)
1645: {
1646: PetscFunctionBegin;
1648: PetscCheck(tao->numbermonitors < MAXTAOMONITORS, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Cannot attach another monitor -- max=%d", MAXTAOMONITORS);
1649: for (PetscInt i = 0; i < tao->numbermonitors; i++) {
1650: PetscBool identical;
1652: PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))(PetscVoidFn *)func, ctx, dest, (PetscErrorCode (*)(void))(PetscVoidFn *)tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i], &identical));
1653: if (identical) PetscFunctionReturn(PETSC_SUCCESS);
1654: }
1655: tao->monitor[tao->numbermonitors] = func;
1656: tao->monitorcontext[tao->numbermonitors] = ctx;
1657: tao->monitordestroy[tao->numbermonitors] = dest;
1658: ++tao->numbermonitors;
1659: PetscFunctionReturn(PETSC_SUCCESS);
1660: }
1662: /*@
1663: TaoMonitorCancel - Clears all the monitor functions for a `Tao` object.
1665: Logically Collective
1667: Input Parameter:
1668: . tao - the `Tao` solver context
1670: Options Database Key:
1671: . -tao_monitor_cancel - cancels all monitors that have been hardwired
1672: into a code by calls to `TaoMonitorSet()`, but does not cancel those
1673: set via the options database
1675: Level: advanced
1677: Note:
1678: There is no way to clear one specific monitor from a `Tao` object.
1680: .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1681: @*/
1682: PetscErrorCode TaoMonitorCancel(Tao tao)
1683: {
1684: PetscInt i;
1686: PetscFunctionBegin;
1688: for (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: PetscInt i;
2654: PetscFunctionBegin;
2656: tao->fc = f;
2657: tao->residual = res;
2658: tao->cnorm = cnorm;
2659: tao->step = steplength;
2660: if (!its) {
2661: tao->cnorm0 = cnorm;
2662: tao->gnorm0 = res;
2663: }
2664: PetscCall(VecLockReadPush(tao->solution));
2665: for (i = 0; i < tao->numbermonitors; i++) PetscCall((*tao->monitor[i])(tao, tao->monitorcontext[i]));
2666: PetscCall(VecLockReadPop(tao->solution));
2667: PetscFunctionReturn(PETSC_SUCCESS);
2668: }
2670: /*@
2671: TaoSetConvergenceHistory - Sets the array used to hold the convergence history.
2673: Logically Collective
2675: Input Parameters:
2676: + tao - the `Tao` solver context
2677: . obj - array to hold objective value history
2678: . resid - array to hold residual history
2679: . cnorm - array to hold constraint violation history
2680: . lits - integer array holds the number of linear iterations for each Tao iteration
2681: . na - size of `obj`, `resid`, and `cnorm`
2682: - reset - `PETSC_TRUE` indicates each new minimization resets the history counter to zero,
2683: else it continues storing new values for new minimizations after the old ones
2685: Level: intermediate
2687: Notes:
2688: If set, `Tao` will fill the given arrays with the indicated
2689: information at each iteration. If 'obj','resid','cnorm','lits' are
2690: *all* `NULL` then space (using size `na`, or 1000 if `na` is `PETSC_DECIDE`) is allocated for the history.
2691: If not all are `NULL`, then only the non-`NULL` information categories
2692: will be stored, the others will be ignored.
2694: Any convergence information after iteration number 'na' will not be stored.
2696: This routine is useful, e.g., when running a code for purposes
2697: of accurate performance monitoring, when no I/O should be done
2698: during the section of code that is being timed.
2700: .seealso: [](ch_tao), `TaoGetConvergenceHistory()`
2701: @*/
2702: PetscErrorCode TaoSetConvergenceHistory(Tao tao, PetscReal obj[], PetscReal resid[], PetscReal cnorm[], PetscInt lits[], PetscInt na, PetscBool reset)
2703: {
2704: PetscFunctionBegin;
2706: if (obj) PetscAssertPointer(obj, 2);
2707: if (resid) PetscAssertPointer(resid, 3);
2708: if (cnorm) PetscAssertPointer(cnorm, 4);
2709: if (lits) PetscAssertPointer(lits, 5);
2711: if (na == PETSC_DECIDE || na == PETSC_CURRENT) na = 1000;
2712: if (!obj && !resid && !cnorm && !lits) {
2713: PetscCall(PetscCalloc4(na, &obj, na, &resid, na, &cnorm, na, &lits));
2714: tao->hist_malloc = PETSC_TRUE;
2715: }
2717: tao->hist_obj = obj;
2718: tao->hist_resid = resid;
2719: tao->hist_cnorm = cnorm;
2720: tao->hist_lits = lits;
2721: tao->hist_max = na;
2722: tao->hist_reset = reset;
2723: tao->hist_len = 0;
2724: PetscFunctionReturn(PETSC_SUCCESS);
2725: }
2727: /*@C
2728: TaoGetConvergenceHistory - Gets the arrays used that hold the convergence history.
2730: Collective
2732: Input Parameter:
2733: . tao - the `Tao` context
2735: Output Parameters:
2736: + obj - array used to hold objective value history
2737: . resid - array used to hold residual history
2738: . cnorm - array used to hold constraint violation history
2739: . lits - integer array used to hold linear solver iteration count
2740: - nhist - size of `obj`, `resid`, `cnorm`, and `lits`
2742: Level: advanced
2744: Notes:
2745: This routine must be preceded by calls to `TaoSetConvergenceHistory()`
2746: and `TaoSolve()`, otherwise it returns useless information.
2748: This routine is useful, e.g., when running a code for purposes
2749: of accurate performance monitoring, when no I/O should be done
2750: during the section of code that is being timed.
2752: Fortran Notes:
2753: The calling sequence is
2754: .vb
2755: call TaoGetConvergenceHistory(Tao tao, PetscInt nhist, PetscErrorCode ierr)
2756: .ve
2757: In other words this gets the current number of entries in the history. Access the history through the array you passed to `TaoSetConvergenceHistory()`
2759: .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetConvergenceHistory()`
2760: @*/
2761: PetscErrorCode TaoGetConvergenceHistory(Tao tao, PetscReal **obj, PetscReal **resid, PetscReal **cnorm, PetscInt **lits, PetscInt *nhist)
2762: {
2763: PetscFunctionBegin;
2765: if (obj) *obj = tao->hist_obj;
2766: if (cnorm) *cnorm = tao->hist_cnorm;
2767: if (resid) *resid = tao->hist_resid;
2768: if (lits) *lits = tao->hist_lits;
2769: if (nhist) *nhist = tao->hist_len;
2770: PetscFunctionReturn(PETSC_SUCCESS);
2771: }
2773: /*@
2774: TaoSetApplicationContext - Sets the optional user-defined context for a `Tao` solver that can be accessed later, for example in the
2775: `Tao` callback functions with `TaoGetApplicationContext()`
2777: Logically Collective
2779: Input Parameters:
2780: + tao - the `Tao` context
2781: - ctx - the user context
2783: Level: intermediate
2785: Fortran Note:
2786: This only works when `ctx` is a Fortran derived type (it cannot be a `PetscObject`), we recommend writing a Fortran interface definition for this
2787: function that tells the Fortran compiler the derived data type that is passed in as the `ctx` argument. See `TaoGetApplicationContext()` for
2788: an example.
2790: .seealso: [](ch_tao), `Tao`, `TaoGetApplicationContext()`
2791: @*/
2792: PetscErrorCode TaoSetApplicationContext(Tao tao, PetscCtx ctx)
2793: {
2794: PetscFunctionBegin;
2796: tao->ctx = ctx;
2797: PetscFunctionReturn(PETSC_SUCCESS);
2798: }
2800: /*@
2801: TaoGetApplicationContext - Gets the user-defined context for a `Tao` solver provided with `TaoSetApplicationContext()`
2803: Not Collective
2805: Input Parameter:
2806: . tao - the `Tao` context
2808: Output Parameter:
2809: . ctx - a pointer to the user context
2811: Level: intermediate
2813: Fortran Note:
2814: This only works when the context is a Fortran derived type or a `PetscObject`. Define `ctx` with
2815: .vb
2816: type(tUsertype), pointer :: ctx
2817: .ve
2819: .seealso: [](ch_tao), `Tao`, `TaoSetApplicationContext()`
2820: @*/
2821: PetscErrorCode TaoGetApplicationContext(Tao tao, PetscCtxRt ctx)
2822: {
2823: PetscFunctionBegin;
2825: PetscAssertPointer(ctx, 2);
2826: *(void **)ctx = tao->ctx;
2827: PetscFunctionReturn(PETSC_SUCCESS);
2828: }
2830: /*@
2831: TaoSetGradientNorm - Sets the matrix used to define the norm that measures the size of the gradient in some of the `Tao` algorithms
2833: Collective
2835: Input Parameters:
2836: + tao - the `Tao` context
2837: - M - matrix that defines the norm
2839: Level: beginner
2841: .seealso: [](ch_tao), `Tao`, `TaoGetGradientNorm()`, `TaoGradientNorm()`
2842: @*/
2843: PetscErrorCode TaoSetGradientNorm(Tao tao, Mat M)
2844: {
2845: PetscFunctionBegin;
2848: PetscCall(PetscObjectReference((PetscObject)M));
2849: PetscCall(MatDestroy(&tao->gradient_norm));
2850: PetscCall(VecDestroy(&tao->gradient_norm_tmp));
2851: tao->gradient_norm = M;
2852: PetscCall(MatCreateVecs(M, NULL, &tao->gradient_norm_tmp));
2853: PetscFunctionReturn(PETSC_SUCCESS);
2854: }
2856: /*@
2857: TaoGetGradientNorm - Returns the matrix used to define the norm used for measuring the size of the gradient in some of the `Tao` algorithms
2859: Not Collective
2861: Input Parameter:
2862: . tao - the `Tao` context
2864: Output Parameter:
2865: . M - gradient norm
2867: Level: beginner
2869: .seealso: [](ch_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGradientNorm()`
2870: @*/
2871: PetscErrorCode TaoGetGradientNorm(Tao tao, Mat *M)
2872: {
2873: PetscFunctionBegin;
2875: PetscAssertPointer(M, 2);
2876: *M = tao->gradient_norm;
2877: PetscFunctionReturn(PETSC_SUCCESS);
2878: }
2880: /*@
2881: TaoGradientNorm - Compute the norm using the `NormType`, the user has selected
2883: Collective
2885: Input Parameters:
2886: + tao - the `Tao` context
2887: . gradient - the gradient
2888: - type - the norm type
2890: Output Parameter:
2891: . gnorm - the gradient norm
2893: Level: advanced
2895: Note:
2896: If `TaoSetGradientNorm()` has been set and `type` is `NORM_2` then the norm provided with `TaoSetGradientNorm()` is used.
2898: Developer Notes:
2899: Should be named `TaoComputeGradientNorm()`.
2901: The usage is a bit confusing, with `TaoSetGradientNorm()` plus `NORM_2` resulting in the computation of the user provided
2902: norm, perhaps a refactorization is in order.
2904: .seealso: [](ch_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGetGradientNorm()`
2905: @*/
2906: PetscErrorCode TaoGradientNorm(Tao tao, Vec gradient, NormType type, PetscReal *gnorm)
2907: {
2908: PetscFunctionBegin;
2912: PetscAssertPointer(gnorm, 4);
2913: if (tao->gradient_norm) {
2914: PetscScalar gnorms;
2916: 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.");
2917: PetscCall(MatMult(tao->gradient_norm, gradient, tao->gradient_norm_tmp));
2918: PetscCall(VecDot(gradient, tao->gradient_norm_tmp, &gnorms));
2919: *gnorm = PetscRealPart(PetscSqrtScalar(gnorms));
2920: } else {
2921: PetscCall(VecNorm(gradient, type, gnorm));
2922: }
2923: PetscFunctionReturn(PETSC_SUCCESS);
2924: }
2926: /*@C
2927: TaoMonitorDrawCtxCreate - Creates the monitor context for `TaoMonitorSolutionDraw()`
2929: Collective
2931: Input Parameters:
2932: + comm - the communicator to share the context
2933: . host - the name of the X Windows host that will display the monitor
2934: . label - the label to put at the top of the display window
2935: . x - the horizontal coordinate of the lower left corner of the window to open
2936: . y - the vertical coordinate of the lower left corner of the window to open
2937: . m - the width of the window
2938: . n - the height of the window
2939: - howoften - how many `Tao` iterations between displaying the monitor information
2941: Output Parameter:
2942: . ctx - the monitor context
2944: Options Database Keys:
2945: + -tao_monitor_solution_draw - use `TaoMonitorSolutionDraw()` to monitor the solution
2946: - -tao_draw_solution_initial - show initial guess as well as current solution
2948: Level: intermediate
2950: Note:
2951: The context this creates, along with `TaoMonitorSolutionDraw()`, and `TaoMonitorDrawCtxDestroy()`
2952: are passed to `TaoMonitorSet()`.
2954: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorDrawCtx()`
2955: @*/
2956: PetscErrorCode TaoMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TaoMonitorDrawCtx *ctx)
2957: {
2958: PetscFunctionBegin;
2959: PetscCall(PetscNew(ctx));
2960: PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
2961: PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
2962: (*ctx)->howoften = howoften;
2963: PetscFunctionReturn(PETSC_SUCCESS);
2964: }
2966: /*@C
2967: TaoMonitorDrawCtxDestroy - Destroys the monitor context for `TaoMonitorSolutionDraw()`
2969: Collective
2971: Input Parameter:
2972: . ictx - the monitor context
2974: Level: intermediate
2976: Note:
2977: This is passed to `TaoMonitorSet()` as the final argument, along with `TaoMonitorSolutionDraw()`, and the context
2978: obtained with `TaoMonitorDrawCtxCreate()`.
2980: .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorSolutionDraw()`
2981: @*/
2982: PetscErrorCode TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx *ictx)
2983: {
2984: PetscFunctionBegin;
2985: PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
2986: PetscCall(PetscFree(*ictx));
2987: PetscFunctionReturn(PETSC_SUCCESS);
2988: }
2990: /*@
2991: TaoGetTerm - Get the entire objective function of the `Tao` as a
2992: single `TaoTerm` in the form $\alpha f(Ax; p)$, where $\alpha$ is a scaling
2993: coefficient, $f$ is a `TaoTerm`, $A$ is an (optional) map and $p$ are the parameters of $f$.
2995: Not collective
2997: Input Parameter:
2998: . tao - a `Tao` context
3000: Output Parameters:
3001: + scale - the scale of the term
3002: . term - a `TaoTerm` for the real-valued function defining the objective
3003: . params - the vector of parameters for `term`, or `NULL` if no parameters were specified for `term`
3004: - map - a map from the solution space of `tao` to the solution space of `term`, if `NULL` then the map is the identity
3006: Level: intermediate
3008: Notes:
3009: If the objective function was defined by providing function callbacks directly to `Tao` (for example, with `TaoSetObjectiveAndGradient()`), then
3010: `TaoGetTerm` will return a `TaoTerm` with the type `TAOTERMCALLBACKS` that encapsulates
3011: those functions.
3013: 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.
3015: .seealso: [](ch_tao), `Tao`, `TaoTerm`, `TAOTERMSUM`, `TaoAddTerm()`
3016: @*/
3017: PetscErrorCode TaoGetTerm(Tao tao, PetscReal *scale, TaoTerm *term, Vec *params, Mat *map)
3018: {
3019: PetscFunctionBegin;
3021: if (scale) PetscAssertPointer(scale, 2);
3022: if (term) PetscAssertPointer(term, 3);
3023: if (params) PetscAssertPointer(params, 4);
3024: if (map) PetscAssertPointer(map, 5);
3025: PetscCall(TaoTermMappingGetData(&tao->objective_term, NULL, scale, term, map));
3026: if (params) *params = tao->objective_parameters;
3027: PetscFunctionReturn(PETSC_SUCCESS);
3028: }
3030: /*@
3031: TaoAddTerm - Add a `term` to the objective function. If `Tao` is empty,
3032: `term` will be the objective of `Tao`.
3034: Collective
3036: Input Parameters:
3037: + tao - a `Tao` solver context
3038: . 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.)
3039: . scale - scaling coefficient for the new term
3040: . term - the real-valued function defining the new term
3041: . params - (optional) parameters for the new term. It is up to each implementation of `TaoTerm` to determine how it behaves when parameters are omitted.
3042: - map - (optional) a map from the `tao` solution space to the `term` solution space; if `NULL` the map is assumed to be the identity
3044: Level: beginner
3046: Notes:
3047: If the objective function was $f(x)$, after calling `TaoAddTerm()` it becomes
3048: $f(x) + \alpha g(Ax; p)$, where $\alpha$ is the `scale`, $g$ is the `term`, $A$ is the
3049: (optional) `map`, and $p$ are the (optional) `params` of $g$.
3051: The `map` $A$ transforms the `Tao` solution vector into the term's solution space.
3052: For example, if the `Tao` solution vector is $x \in \mathbb{R}^n$ and the mapping
3053: matrix is $A \in \mathbb{R}^{m \times n}$, then the term evaluates $g(Ax; p)$ with
3054: $Ax \in \mathbb{R}^m$. The term's solution space is therefore $\mathbb{R}^m$. If the map is
3055: `NULL`, the identity is used and the term's solution space must match the `Tao` solution space.
3056: `Tao` automatically applies the chain rule for gradients ($A^T \nabla g$) and Hessians
3057: ($A^T \nabla^2 g \, A$) with respect to $x$.
3059: The `params` $p$ are fixed data that are not optimized over. Some `TaoTermType`s
3060: require the parameter space to be related to the term's solution space (e.g., the same
3061: size); when a mapping matrix $A$ is used, the parameter space may depend on either the row
3062: or column space of $A$. See the documentation for each `TaoTermType`.
3064: Currently, `TaoAddTerm()` does not support bounded Newton solvers (`TAOBNK`,`TAOBNLS`,`TAOBNTL`,`TAOBNTR`,and `TAOBQNK`)
3066: .seealso: [](ch_tao), `Tao`, `TaoTerm`, `TAOTERMSUM`, `TaoGetTerm()`
3067: @*/
3068: PetscErrorCode TaoAddTerm(Tao tao, const char prefix[], PetscReal scale, TaoTerm term, Vec params, Mat map)
3069: {
3070: PetscBool is_sum, is_callback;
3071: PetscInt num_old_terms;
3072: Vec *vec_list = NULL;
3074: PetscFunctionBegin;
3076: if (prefix) PetscAssertPointer(prefix, 2);
3079: PetscCheckSameComm(tao, 1, term, 4);
3080: if (params) {
3082: PetscCheckSameComm(tao, 1, params, 5);
3083: }
3084: if (map) {
3086: PetscCheckSameComm(tao, 1, map, 6);
3087: }
3088: // If user is using TaoAddTerm, before setting any terms or callbacks,
3089: // then tao->objective_term.term is empty callback, which we want to remove.
3090: PetscCall(PetscObjectTypeCompare((PetscObject)tao->objective_term.term, TAOTERMCALLBACKS, &is_callback));
3091: PetscCall(PetscObjectTypeCompare((PetscObject)term, TAOTERMSUM, &is_sum));
3092: PetscCheck(!is_sum, PetscObjectComm((PetscObject)term), PETSC_ERR_ARG_WRONG, "TaoAddTerm does not support adding TAOTERMSUM");
3093: if (is_callback) {
3094: PetscBool is_obj, is_objgrad, is_grad;
3096: PetscCall(TaoTermIsObjectiveDefined(tao->objective_term.term, &is_obj));
3097: PetscCall(TaoTermIsObjectiveAndGradientDefined(tao->objective_term.term, &is_objgrad));
3098: PetscCall(TaoTermIsGradientDefined(tao->objective_term.term, &is_grad));
3099: // Empty callback term
3100: if (!(is_obj || is_objgrad || is_grad)) {
3101: PetscCall(TaoTermMappingSetData(&tao->objective_term, NULL, scale, term, map));
3102: PetscCall(PetscObjectReference((PetscObject)params));
3103: PetscCall(VecDestroy(&tao->objective_parameters));
3104: // Empty callback term. Destroy hessians, as they are not needed
3105: PetscCall(MatDestroy(&tao->hessian));
3106: PetscCall(MatDestroy(&tao->hessian_pre));
3107: tao->objective_parameters = params;
3108: tao->term_set = PETSC_TRUE;
3109: PetscFunctionReturn(PETSC_SUCCESS);
3110: }
3111: }
3112: PetscCall(PetscObjectTypeCompare((PetscObject)tao->objective_term.term, TAOTERMSUM, &is_sum));
3113: // One TaoTerm has been set. Create TAOTERMSUM to store that, and the new one
3114: if (!is_sum) {
3115: TaoTerm old_sum;
3116: const char *tao_prefix;
3117: const char *term_prefix;
3119: PetscCall(TaoTermDuplicate(tao->objective_term.term, TAOTERM_DUPLICATE_SIZEONLY, &old_sum));
3120: if (tao->objective_term.map) {
3121: VecType map_vectype;
3122: VecType param_vectype;
3123: PetscLayout cmap, param_layout;
3125: PetscCall(MatGetVecType(tao->objective_term.map, &map_vectype));
3126: PetscCall(MatGetLayouts(tao->objective_term.map, NULL, &cmap));
3127: PetscCall(TaoTermGetParametersVecType(old_sum, ¶m_vectype));
3128: PetscCall(TaoTermGetParametersLayout(old_sum, ¶m_layout));
3130: PetscCall(TaoTermSetSolutionVecType(old_sum, map_vectype));
3131: PetscCall(TaoTermSetParametersVecType(old_sum, param_vectype));
3132: PetscCall(TaoTermSetSolutionLayout(old_sum, cmap));
3133: PetscCall(TaoTermSetParametersLayout(old_sum, param_layout));
3134: }
3136: PetscCall(TaoTermSetType(old_sum, TAOTERMSUM));
3137: PetscCall(TaoGetOptionsPrefix(tao, &tao_prefix));
3138: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)old_sum, tao_prefix));
3139: PetscCall(TaoTermSumSetNumberTerms(old_sum, 1));
3140: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)tao->objective_term.term, &term_prefix));
3141: PetscCall(TaoTermSumSetTerm(old_sum, 0, term_prefix, tao->objective_term.scale, tao->objective_term.term, tao->objective_term.map));
3142: PetscCall(TaoTermSumSetTermHessianMatrices(old_sum, 0, NULL, NULL, tao->hessian, tao->hessian_pre));
3143: PetscCall(MatDestroy(&tao->hessian));
3144: PetscCall(MatDestroy(&tao->hessian_pre));
3145: PetscCall(TaoTermMappingReset(&tao->objective_term));
3146: PetscCall(TaoTermMappingSetData(&tao->objective_term, NULL, 1.0, old_sum, NULL));
3147: if (tao->objective_parameters) {
3148: // convert the parameters to a VECNEST
3149: Vec subvecs[1];
3151: subvecs[0] = tao->objective_parameters;
3152: tao->objective_parameters = NULL;
3153: PetscCall(TaoTermSumParametersPack(old_sum, subvecs, &tao->objective_parameters));
3154: PetscCall(VecDestroy(&subvecs[0]));
3155: }
3156: PetscCall(TaoTermDestroy(&old_sum));
3157: tao->num_terms = 1;
3158: }
3159: PetscCall(TaoTermSumGetNumberTerms(tao->objective_term.term, &num_old_terms));
3160: if (tao->objective_parameters || params) {
3161: PetscCall(PetscCalloc1(num_old_terms + 1, &vec_list));
3162: if (tao->objective_parameters) PetscCall(TaoTermSumParametersUnpack(tao->objective_term.term, &tao->objective_parameters, vec_list));
3163: PetscCall(PetscObjectReference((PetscObject)params));
3164: vec_list[num_old_terms] = params;
3165: }
3166: PetscCall(TaoTermSumAddTerm(tao->objective_term.term, prefix, scale, term, map, NULL));
3167: tao->num_terms++;
3168: if (vec_list) {
3169: PetscInt num_terms = num_old_terms + 1;
3170: PetscCall(TaoTermSumParametersPack(tao->objective_term.term, vec_list, &tao->objective_parameters));
3171: for (PetscInt i = 0; i < num_terms; i++) PetscCall(VecDestroy(&vec_list[i]));
3172: PetscCall(PetscFree(vec_list));
3173: }
3174: PetscFunctionReturn(PETSC_SUCCESS);
3175: }