1 #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/ 2 #include <petsc/private/snesimpl.h> 3 4 PetscBool TaoRegisterAllCalled = PETSC_FALSE; 5 PetscFunctionList TaoList = NULL; 6 7 PetscClassId TAO_CLASSID; 8 9 PetscLogEvent TAO_Solve; 10 PetscLogEvent TAO_ObjectiveEval; 11 PetscLogEvent TAO_GradientEval; 12 PetscLogEvent TAO_ObjGradEval; 13 PetscLogEvent TAO_HessianEval; 14 PetscLogEvent TAO_JacobianEval; 15 PetscLogEvent TAO_ConstraintsEval; 16 17 const char *TaoSubSetTypes[] = {"subvec", "mask", "matrixfree", "TaoSubSetType", "TAO_SUBSET_", NULL}; 18 19 struct _n_TaoMonitorDrawCtx { 20 PetscViewer viewer; 21 PetscInt howoften; /* when > 0 uses iteration % howoften, when negative only final solution plotted */ 22 }; 23 24 static PetscErrorCode KSPPreSolve_TAOEW_Private(KSP ksp, Vec b, Vec x, Tao tao) 25 { 26 SNES snes_ewdummy = tao->snes_ewdummy; 27 28 PetscFunctionBegin; 29 if (!snes_ewdummy) PetscFunctionReturn(PETSC_SUCCESS); 30 /* populate snes_ewdummy struct values used in KSPPreSolve_SNESEW */ 31 snes_ewdummy->vec_func = b; 32 snes_ewdummy->rtol = tao->gttol; 33 snes_ewdummy->iter = tao->niter; 34 PetscCall(VecNorm(b, NORM_2, &snes_ewdummy->norm)); 35 PetscCall(KSPPreSolve_SNESEW(ksp, b, x, snes_ewdummy)); 36 snes_ewdummy->vec_func = NULL; 37 PetscFunctionReturn(PETSC_SUCCESS); 38 } 39 40 static PetscErrorCode KSPPostSolve_TAOEW_Private(KSP ksp, Vec b, Vec x, Tao tao) 41 { 42 SNES snes_ewdummy = tao->snes_ewdummy; 43 44 PetscFunctionBegin; 45 if (!snes_ewdummy) PetscFunctionReturn(PETSC_SUCCESS); 46 PetscCall(KSPPostSolve_SNESEW(ksp, b, x, snes_ewdummy)); 47 PetscFunctionReturn(PETSC_SUCCESS); 48 } 49 50 static PetscErrorCode TaoSetUpEW_Private(Tao tao) 51 { 52 SNESKSPEW *kctx; 53 const char *ewprefix; 54 55 PetscFunctionBegin; 56 if (!tao->ksp) PetscFunctionReturn(PETSC_SUCCESS); 57 if (tao->ksp_ewconv) { 58 if (!tao->snes_ewdummy) PetscCall(SNESCreate(PetscObjectComm((PetscObject)tao), &tao->snes_ewdummy)); 59 tao->snes_ewdummy->ksp_ewconv = PETSC_TRUE; 60 PetscCall(KSPSetPreSolve(tao->ksp, (PetscErrorCode(*)(KSP, Vec, Vec, void *))KSPPreSolve_TAOEW_Private, tao)); 61 PetscCall(KSPSetPostSolve(tao->ksp, (PetscErrorCode(*)(KSP, Vec, Vec, void *))KSPPostSolve_TAOEW_Private, tao)); 62 63 PetscCall(KSPGetOptionsPrefix(tao->ksp, &ewprefix)); 64 kctx = (SNESKSPEW *)tao->snes_ewdummy->kspconvctx; 65 PetscCall(SNESEWSetFromOptions_Private(kctx, PETSC_FALSE, PetscObjectComm((PetscObject)tao), ewprefix)); 66 } else PetscCall(SNESDestroy(&tao->snes_ewdummy)); 67 PetscFunctionReturn(PETSC_SUCCESS); 68 } 69 70 /*@ 71 TaoCreate - Creates a Tao solver 72 73 Collective 74 75 Input Parameter: 76 . comm - MPI communicator 77 78 Output Parameter: 79 . newtao - the new `Tao` context 80 81 Options Database Key: 82 . -tao_type - select which method Tao should use 83 84 Level: beginner 85 86 .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoDestroy()`, `TAOSetFromOptions()`, `TAOSetType()` 87 @*/ 88 PetscErrorCode TaoCreate(MPI_Comm comm, Tao *newtao) 89 { 90 Tao tao; 91 92 PetscFunctionBegin; 93 PetscAssertPointer(newtao, 2); 94 PetscCall(TaoInitializePackage()); 95 PetscCall(TaoLineSearchInitializePackage()); 96 PetscCall(PetscHeaderCreate(tao, TAO_CLASSID, "Tao", "Optimization solver", "Tao", comm, TaoDestroy, TaoView)); 97 98 /* Set non-NULL defaults */ 99 tao->ops->convergencetest = TaoDefaultConvergenceTest; 100 101 tao->max_it = 10000; 102 tao->max_funcs = -1; 103 #if defined(PETSC_USE_REAL_SINGLE) 104 tao->gatol = 1e-5; 105 tao->grtol = 1e-5; 106 tao->crtol = 1e-5; 107 tao->catol = 1e-5; 108 #else 109 tao->gatol = 1e-8; 110 tao->grtol = 1e-8; 111 tao->crtol = 1e-8; 112 tao->catol = 1e-8; 113 #endif 114 tao->gttol = 0.0; 115 tao->steptol = 0.0; 116 tao->trust0 = PETSC_INFINITY; 117 tao->fmin = PETSC_NINFINITY; 118 119 tao->hist_reset = PETSC_TRUE; 120 121 PetscCall(TaoResetStatistics(tao)); 122 *newtao = tao; 123 PetscFunctionReturn(PETSC_SUCCESS); 124 } 125 126 /*@ 127 TaoSolve - Solves an optimization problem min F(x) s.t. l <= x <= u 128 129 Collective 130 131 Input Parameter: 132 . tao - the `Tao` context 133 134 Level: beginner 135 136 Notes: 137 The user must set up the `Tao` object with calls to `TaoSetSolution()`, `TaoSetObjective()`, `TaoSetGradient()`, and (if using 2nd order method) `TaoSetHessian()`. 138 139 You should call `TaoGetConvergedReason()` or run with `-tao_converged_reason` to determine if the optimization algorithm actually succeeded or 140 why it failed. 141 142 .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSetObjective()`, `TaoSetGradient()`, `TaoSetHessian()`, `TaoGetConvergedReason()`, `TaoSetUp()` 143 @*/ 144 PetscErrorCode TaoSolve(Tao tao) 145 { 146 static PetscBool set = PETSC_FALSE; 147 148 PetscFunctionBegin; 149 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 150 PetscCall(PetscCitationsRegister("@TechReport{tao-user-ref,\n" 151 "title = {Toolkit for Advanced Optimization (TAO) Users Manual},\n" 152 "author = {Todd Munson and Jason Sarich and Stefan Wild and Steve Benson and Lois Curfman McInnes},\n" 153 "Institution = {Argonne National Laboratory},\n" 154 "Year = 2014,\n" 155 "Number = {ANL/MCS-TM-322 - Revision 3.5},\n" 156 "url = {https://www.mcs.anl.gov/research/projects/tao/}\n}\n", 157 &set)); 158 tao->header_printed = PETSC_FALSE; 159 PetscCall(TaoSetUp(tao)); 160 PetscCall(TaoResetStatistics(tao)); 161 if (tao->linesearch) PetscCall(TaoLineSearchReset(tao->linesearch)); 162 163 PetscCall(PetscLogEventBegin(TAO_Solve, tao, 0, 0, 0)); 164 PetscTryTypeMethod(tao, solve); 165 PetscCall(PetscLogEventEnd(TAO_Solve, tao, 0, 0, 0)); 166 167 PetscCall(VecViewFromOptions(tao->solution, (PetscObject)tao, "-tao_view_solution")); 168 169 tao->ntotalits += tao->niter; 170 171 if (tao->printreason) { 172 PetscViewer viewer = PETSC_VIEWER_STDOUT_(((PetscObject)tao)->comm); 173 PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)tao)->tablevel)); 174 if (tao->reason > 0) { 175 PetscCall(PetscViewerASCIIPrintf(viewer, " TAO %s solve converged due to %s iterations %" PetscInt_FMT "\n", ((PetscObject)tao)->prefix ? ((PetscObject)tao)->prefix : "", TaoConvergedReasons[tao->reason], tao->niter)); 176 } else { 177 PetscCall(PetscViewerASCIIPrintf(viewer, " TAO %s solve did not converge due to %s iteration %" PetscInt_FMT "\n", ((PetscObject)tao)->prefix ? ((PetscObject)tao)->prefix : "", TaoConvergedReasons[tao->reason], tao->niter)); 178 } 179 PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)tao)->tablevel)); 180 } 181 PetscCall(TaoViewFromOptions(tao, NULL, "-tao_view")); 182 PetscFunctionReturn(PETSC_SUCCESS); 183 } 184 185 /*@ 186 TaoSetUp - Sets up the internal data structures for the later use 187 of a Tao solver 188 189 Collective 190 191 Input Parameter: 192 . tao - the `Tao` context 193 194 Level: advanced 195 196 Note: 197 The user will not need to explicitly call `TaoSetUp()`, as it will 198 automatically be called in `TaoSolve()`. However, if the user 199 desires to call it explicitly, it should come after `TaoCreate()` 200 and any TaoSetSomething() routines, but before `TaoSolve()`. 201 202 .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()` 203 @*/ 204 PetscErrorCode TaoSetUp(Tao tao) 205 { 206 PetscFunctionBegin; 207 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 208 if (tao->setupcalled) PetscFunctionReturn(PETSC_SUCCESS); 209 PetscCall(TaoSetUpEW_Private(tao)); 210 PetscCheck(tao->solution, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoSetSolution"); 211 PetscTryTypeMethod(tao, setup); 212 tao->setupcalled = PETSC_TRUE; 213 PetscFunctionReturn(PETSC_SUCCESS); 214 } 215 216 /*@C 217 TaoDestroy - Destroys the `Tao` context that was created with `TaoCreate()` 218 219 Collective 220 221 Input Parameter: 222 . tao - the `Tao` context 223 224 Level: beginner 225 226 .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()` 227 @*/ 228 PetscErrorCode TaoDestroy(Tao *tao) 229 { 230 PetscFunctionBegin; 231 if (!*tao) PetscFunctionReturn(PETSC_SUCCESS); 232 PetscValidHeaderSpecific(*tao, TAO_CLASSID, 1); 233 if (--((PetscObject)*tao)->refct > 0) { 234 *tao = NULL; 235 PetscFunctionReturn(PETSC_SUCCESS); 236 } 237 238 if ((*tao)->ops->destroy) PetscCall((*((*tao))->ops->destroy)(*tao)); 239 PetscCall(KSPDestroy(&(*tao)->ksp)); 240 PetscCall(SNESDestroy(&(*tao)->snes_ewdummy)); 241 PetscCall(TaoLineSearchDestroy(&(*tao)->linesearch)); 242 243 if ((*tao)->ops->convergencedestroy) { 244 PetscCall((*(*tao)->ops->convergencedestroy)((*tao)->cnvP)); 245 if ((*tao)->jacobian_state_inv) PetscCall(MatDestroy(&(*tao)->jacobian_state_inv)); 246 } 247 PetscCall(VecDestroy(&(*tao)->solution)); 248 PetscCall(VecDestroy(&(*tao)->gradient)); 249 PetscCall(VecDestroy(&(*tao)->ls_res)); 250 251 if ((*tao)->gradient_norm) { 252 PetscCall(PetscObjectDereference((PetscObject)(*tao)->gradient_norm)); 253 PetscCall(VecDestroy(&(*tao)->gradient_norm_tmp)); 254 } 255 256 PetscCall(VecDestroy(&(*tao)->XL)); 257 PetscCall(VecDestroy(&(*tao)->XU)); 258 PetscCall(VecDestroy(&(*tao)->IL)); 259 PetscCall(VecDestroy(&(*tao)->IU)); 260 PetscCall(VecDestroy(&(*tao)->DE)); 261 PetscCall(VecDestroy(&(*tao)->DI)); 262 PetscCall(VecDestroy(&(*tao)->constraints)); 263 PetscCall(VecDestroy(&(*tao)->constraints_equality)); 264 PetscCall(VecDestroy(&(*tao)->constraints_inequality)); 265 PetscCall(VecDestroy(&(*tao)->stepdirection)); 266 PetscCall(MatDestroy(&(*tao)->hessian_pre)); 267 PetscCall(MatDestroy(&(*tao)->hessian)); 268 PetscCall(MatDestroy(&(*tao)->ls_jac)); 269 PetscCall(MatDestroy(&(*tao)->ls_jac_pre)); 270 PetscCall(MatDestroy(&(*tao)->jacobian_pre)); 271 PetscCall(MatDestroy(&(*tao)->jacobian)); 272 PetscCall(MatDestroy(&(*tao)->jacobian_state_pre)); 273 PetscCall(MatDestroy(&(*tao)->jacobian_state)); 274 PetscCall(MatDestroy(&(*tao)->jacobian_state_inv)); 275 PetscCall(MatDestroy(&(*tao)->jacobian_design)); 276 PetscCall(MatDestroy(&(*tao)->jacobian_equality)); 277 PetscCall(MatDestroy(&(*tao)->jacobian_equality_pre)); 278 PetscCall(MatDestroy(&(*tao)->jacobian_inequality)); 279 PetscCall(MatDestroy(&(*tao)->jacobian_inequality_pre)); 280 PetscCall(ISDestroy(&(*tao)->state_is)); 281 PetscCall(ISDestroy(&(*tao)->design_is)); 282 PetscCall(VecDestroy(&(*tao)->res_weights_v)); 283 PetscCall(TaoCancelMonitors(*tao)); 284 if ((*tao)->hist_malloc) PetscCall(PetscFree4((*tao)->hist_obj, (*tao)->hist_resid, (*tao)->hist_cnorm, (*tao)->hist_lits)); 285 if ((*tao)->res_weights_n) { 286 PetscCall(PetscFree((*tao)->res_weights_rows)); 287 PetscCall(PetscFree((*tao)->res_weights_cols)); 288 PetscCall(PetscFree((*tao)->res_weights_w)); 289 } 290 PetscCall(PetscHeaderDestroy(tao)); 291 PetscFunctionReturn(PETSC_SUCCESS); 292 } 293 294 /*@ 295 TaoKSPSetUseEW - Sets `SNES` to use Eisenstat-Walker method {cite}`ew96`for computing relative tolerance for linear solvers. 296 297 Logically Collective 298 299 Input Parameters: 300 + tao - Tao context 301 - flag - `PETSC_TRUE` or `PETSC_FALSE` 302 303 Level: advanced 304 305 Note: 306 See `SNESKSPSetUseEW()` for customization details. 307 308 .seealso: [](ch_tao), `Tao`, `SNESKSPSetUseEW()` 309 @*/ 310 PetscErrorCode TaoKSPSetUseEW(Tao tao, PetscBool flag) 311 { 312 PetscFunctionBegin; 313 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 314 PetscValidLogicalCollectiveBool(tao, flag, 2); 315 tao->ksp_ewconv = flag; 316 PetscFunctionReturn(PETSC_SUCCESS); 317 } 318 319 /*@ 320 TaoSetFromOptions - Sets various Tao parameters from the options database 321 322 Collective 323 324 Input Parameter: 325 . tao - the `Tao` solver context 326 327 Options Database Keys: 328 + -tao_type <type> - The algorithm that Tao uses (lmvm, nls, etc.) 329 . -tao_gatol <gatol> - absolute error tolerance for ||gradient|| 330 . -tao_grtol <grtol> - relative error tolerance for ||gradient|| 331 . -tao_gttol <gttol> - reduction of ||gradient|| relative to initial gradient 332 . -tao_max_it <max> - sets maximum number of iterations 333 . -tao_max_funcs <max> - sets maximum number of function evaluations 334 . -tao_fmin <fmin> - stop if function value reaches fmin 335 . -tao_steptol <tol> - stop if trust region radius less than <tol> 336 . -tao_trust0 <t> - initial trust region radius 337 . -tao_monitor - prints function value and residual norm at each iteration 338 . -tao_smonitor - same as tao_monitor, but truncates very small values 339 . -tao_cmonitor - prints function value, residual, and constraint norm at each iteration 340 . -tao_view_solution - prints solution vector at each iteration 341 . -tao_view_ls_residual - prints least-squares residual vector at each iteration 342 . -tao_view_stepdirection - prints step direction vector at each iteration 343 . -tao_view_gradient - prints gradient vector at each iteration 344 . -tao_draw_solution - graphically view solution vector at each iteration 345 . -tao_draw_step - graphically view step vector at each iteration 346 . -tao_draw_gradient - graphically view gradient at each iteration 347 . -tao_fd_gradient - use gradient computed with finite differences 348 . -tao_fd_hessian - use hessian computed with finite differences 349 . -tao_mf_hessian - use matrix-free hessian computed with finite differences 350 . -tao_cancelmonitors - cancels all monitors (except those set with command line) 351 . -tao_view - prints information about the Tao after solving 352 - -tao_converged_reason - prints the reason Tao stopped iterating 353 354 Level: beginner 355 356 Note: 357 To see all options, run your program with the `-help` option or consult the 358 user's manual. Should be called after `TaoCreate()` but before `TaoSolve()` 359 360 .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()` 361 @*/ 362 PetscErrorCode TaoSetFromOptions(Tao tao) 363 { 364 TaoType default_type = TAOLMVM; 365 char type[256], monfilename[PETSC_MAX_PATH_LEN]; 366 PetscViewer monviewer; 367 PetscBool flg; 368 MPI_Comm comm; 369 370 PetscFunctionBegin; 371 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 372 PetscCall(PetscObjectGetComm((PetscObject)tao, &comm)); 373 374 if (((PetscObject)tao)->type_name) default_type = ((PetscObject)tao)->type_name; 375 376 PetscObjectOptionsBegin((PetscObject)tao); 377 /* Check for type from options */ 378 PetscCall(PetscOptionsFList("-tao_type", "Tao Solver type", "TaoSetType", TaoList, default_type, type, 256, &flg)); 379 if (flg) { 380 PetscCall(TaoSetType(tao, type)); 381 } else if (!((PetscObject)tao)->type_name) { 382 PetscCall(TaoSetType(tao, default_type)); 383 } 384 385 /* Tao solvers do not set the prefix, set it here if not yet done 386 We do it after SetType since solver may have been changed */ 387 if (tao->linesearch) { 388 const char *prefix; 389 PetscCall(TaoLineSearchGetOptionsPrefix(tao->linesearch, &prefix)); 390 if (!prefix) PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, ((PetscObject)(tao))->prefix)); 391 } 392 393 PetscCall(PetscOptionsReal("-tao_catol", "Stop if constraints violations within", "TaoSetConstraintTolerances", tao->catol, &tao->catol, &flg)); 394 if (flg) tao->catol_changed = PETSC_TRUE; 395 PetscCall(PetscOptionsReal("-tao_crtol", "Stop if relative constraint violations within", "TaoSetConstraintTolerances", tao->crtol, &tao->crtol, &flg)); 396 if (flg) tao->crtol_changed = PETSC_TRUE; 397 PetscCall(PetscOptionsReal("-tao_gatol", "Stop if norm of gradient less than", "TaoSetTolerances", tao->gatol, &tao->gatol, &flg)); 398 if (flg) tao->gatol_changed = PETSC_TRUE; 399 PetscCall(PetscOptionsReal("-tao_grtol", "Stop if norm of gradient divided by the function value is less than", "TaoSetTolerances", tao->grtol, &tao->grtol, &flg)); 400 if (flg) tao->grtol_changed = PETSC_TRUE; 401 PetscCall(PetscOptionsReal("-tao_gttol", "Stop if the norm of the gradient is less than the norm of the initial gradient times tol", "TaoSetTolerances", tao->gttol, &tao->gttol, &flg)); 402 if (flg) tao->gttol_changed = PETSC_TRUE; 403 PetscCall(PetscOptionsInt("-tao_max_it", "Stop if iteration number exceeds", "TaoSetMaximumIterations", tao->max_it, &tao->max_it, &flg)); 404 if (flg) tao->max_it_changed = PETSC_TRUE; 405 PetscCall(PetscOptionsInt("-tao_max_funcs", "Stop if number of function evaluations exceeds", "TaoSetMaximumFunctionEvaluations", tao->max_funcs, &tao->max_funcs, &flg)); 406 if (flg) tao->max_funcs_changed = PETSC_TRUE; 407 PetscCall(PetscOptionsReal("-tao_fmin", "Stop if function less than", "TaoSetFunctionLowerBound", tao->fmin, &tao->fmin, &flg)); 408 if (flg) tao->fmin_changed = PETSC_TRUE; 409 PetscCall(PetscOptionsReal("-tao_steptol", "Stop if step size or trust region radius less than", "", tao->steptol, &tao->steptol, &flg)); 410 if (flg) tao->steptol_changed = PETSC_TRUE; 411 PetscCall(PetscOptionsReal("-tao_trust0", "Initial trust region radius", "TaoSetTrustRegionRadius", tao->trust0, &tao->trust0, &flg)); 412 if (flg) tao->trust0_changed = PETSC_TRUE; 413 PetscCall(PetscOptionsString("-tao_view_solution", "view solution vector after each evaluation", "TaoSetMonitor", "stdout", monfilename, sizeof(monfilename), &flg)); 414 if (flg) { 415 PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer)); 416 PetscCall(TaoSetMonitor(tao, TaoSolutionMonitor, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy)); 417 } 418 419 PetscCall(PetscOptionsBool("-tao_converged_reason", "Print reason for Tao converged", "TaoSolve", tao->printreason, &tao->printreason, NULL)); 420 PetscCall(PetscOptionsString("-tao_view_gradient", "view gradient vector after each evaluation", "TaoSetMonitor", "stdout", monfilename, sizeof(monfilename), &flg)); 421 if (flg) { 422 PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer)); 423 PetscCall(TaoSetMonitor(tao, TaoGradientMonitor, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy)); 424 } 425 426 PetscCall(PetscOptionsString("-tao_view_stepdirection", "view step direction vector after each iteration", "TaoSetMonitor", "stdout", monfilename, sizeof(monfilename), &flg)); 427 if (flg) { 428 PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer)); 429 PetscCall(TaoSetMonitor(tao, TaoStepDirectionMonitor, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy)); 430 } 431 432 PetscCall(PetscOptionsString("-tao_view_residual", "view least-squares residual vector after each evaluation", "TaoSetMonitor", "stdout", monfilename, sizeof(monfilename), &flg)); 433 if (flg) { 434 PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer)); 435 PetscCall(TaoSetMonitor(tao, TaoResidualMonitor, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy)); 436 } 437 438 PetscCall(PetscOptionsString("-tao_monitor", "Use the default convergence monitor", "TaoSetMonitor", "stdout", monfilename, sizeof(monfilename), &flg)); 439 if (flg) { 440 PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer)); 441 PetscCall(TaoSetMonitor(tao, TaoMonitorDefault, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy)); 442 } 443 444 PetscCall(PetscOptionsString("-tao_gmonitor", "Use the convergence monitor with extra globalization info", "TaoSetMonitor", "stdout", monfilename, sizeof(monfilename), &flg)); 445 if (flg) { 446 PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer)); 447 PetscCall(TaoSetMonitor(tao, TaoDefaultGMonitor, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy)); 448 } 449 450 PetscCall(PetscOptionsString("-tao_smonitor", "Use the short convergence monitor", "TaoSetMonitor", "stdout", monfilename, sizeof(monfilename), &flg)); 451 if (flg) { 452 PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer)); 453 PetscCall(TaoSetMonitor(tao, TaoDefaultSMonitor, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy)); 454 } 455 456 PetscCall(PetscOptionsString("-tao_cmonitor", "Use the default convergence monitor with constraint norm", "TaoSetMonitor", "stdout", monfilename, sizeof(monfilename), &flg)); 457 if (flg) { 458 PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer)); 459 PetscCall(TaoSetMonitor(tao, TaoDefaultCMonitor, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy)); 460 } 461 462 flg = PETSC_FALSE; 463 PetscCall(PetscOptionsBool("-tao_cancelmonitors", "cancel all monitors and call any registered destroy routines", "TaoCancelMonitors", flg, &flg, NULL)); 464 if (flg) PetscCall(TaoCancelMonitors(tao)); 465 466 flg = PETSC_FALSE; 467 PetscCall(PetscOptionsBool("-tao_draw_solution", "Plot solution vector at each iteration", "TaoSetMonitor", flg, &flg, NULL)); 468 if (flg) { 469 TaoMonitorDrawCtx drawctx; 470 PetscInt howoften = 1; 471 PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &drawctx)); 472 PetscCall(TaoSetMonitor(tao, TaoDrawSolutionMonitor, drawctx, (PetscErrorCode(*)(void **))TaoMonitorDrawCtxDestroy)); 473 } 474 475 flg = PETSC_FALSE; 476 PetscCall(PetscOptionsBool("-tao_draw_step", "plots step direction at each iteration", "TaoSetMonitor", flg, &flg, NULL)); 477 if (flg) PetscCall(TaoSetMonitor(tao, TaoDrawStepMonitor, NULL, NULL)); 478 479 flg = PETSC_FALSE; 480 PetscCall(PetscOptionsBool("-tao_draw_gradient", "plots gradient at each iteration", "TaoSetMonitor", flg, &flg, NULL)); 481 if (flg) { 482 TaoMonitorDrawCtx drawctx; 483 PetscInt howoften = 1; 484 PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &drawctx)); 485 PetscCall(TaoSetMonitor(tao, TaoDrawGradientMonitor, drawctx, (PetscErrorCode(*)(void **))TaoMonitorDrawCtxDestroy)); 486 } 487 flg = PETSC_FALSE; 488 PetscCall(PetscOptionsBool("-tao_fd_gradient", "compute gradient using finite differences", "TaoDefaultComputeGradient", flg, &flg, NULL)); 489 if (flg) PetscCall(TaoSetGradient(tao, NULL, TaoDefaultComputeGradient, NULL)); 490 flg = PETSC_FALSE; 491 PetscCall(PetscOptionsBool("-tao_fd_hessian", "compute hessian using finite differences", "TaoDefaultComputeHessian", flg, &flg, NULL)); 492 if (flg) { 493 Mat H; 494 495 PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &H)); 496 PetscCall(MatSetType(H, MATAIJ)); 497 PetscCall(TaoSetHessian(tao, H, H, TaoDefaultComputeHessian, NULL)); 498 PetscCall(MatDestroy(&H)); 499 } 500 flg = PETSC_FALSE; 501 PetscCall(PetscOptionsBool("-tao_mf_hessian", "compute matrix-free hessian using finite differences", "TaoDefaultComputeHessianMFFD", flg, &flg, NULL)); 502 if (flg) { 503 Mat H; 504 505 PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &H)); 506 PetscCall(TaoSetHessian(tao, H, H, TaoDefaultComputeHessianMFFD, NULL)); 507 PetscCall(MatDestroy(&H)); 508 } 509 flg = PETSC_FALSE; 510 PetscCall(PetscOptionsBool("-tao_recycle_history", "enable recycling/re-using information from the previous TaoSolve() call for some algorithms", "TaoSetRecycleHistory", flg, &flg, NULL)); 511 if (flg) PetscCall(TaoSetRecycleHistory(tao, PETSC_TRUE)); 512 PetscCall(PetscOptionsEnum("-tao_subset_type", "subset type", "", TaoSubSetTypes, (PetscEnum)tao->subset_type, (PetscEnum *)&tao->subset_type, NULL)); 513 514 if (tao->ksp) { 515 PetscCall(PetscOptionsBool("-tao_ksp_ew", "Use Eisentat-Walker linear system convergence test", "TaoKSPSetUseEW", tao->ksp_ewconv, &tao->ksp_ewconv, NULL)); 516 PetscCall(TaoKSPSetUseEW(tao, tao->ksp_ewconv)); 517 } 518 519 PetscTryTypeMethod(tao, setfromoptions, PetscOptionsObject); 520 521 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 522 PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)tao, PetscOptionsObject)); 523 PetscOptionsEnd(); 524 525 if (tao->linesearch) PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 526 PetscFunctionReturn(PETSC_SUCCESS); 527 } 528 529 /*@C 530 TaoViewFromOptions - View a `Tao` object based on values in the options database 531 532 Collective 533 534 Input Parameters: 535 + A - the `Tao` context 536 . obj - Optional object that provides the prefix for the options database 537 - name - command line option 538 539 Level: intermediate 540 541 .seealso: [](ch_tao), `Tao`, `TaoView`, `PetscObjectViewFromOptions()`, `TaoCreate()` 542 @*/ 543 PetscErrorCode TaoViewFromOptions(Tao A, PetscObject obj, const char name[]) 544 { 545 PetscFunctionBegin; 546 PetscValidHeaderSpecific(A, TAO_CLASSID, 1); 547 PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 548 PetscFunctionReturn(PETSC_SUCCESS); 549 } 550 551 /*@C 552 TaoView - Prints information about the `Tao` object 553 554 Collective 555 556 Input Parameters: 557 + tao - the `Tao` context 558 - viewer - visualization context 559 560 Options Database Key: 561 . -tao_view - Calls `TaoView()` at the end of `TaoSolve()` 562 563 Level: beginner 564 565 Notes: 566 The available visualization contexts include 567 + `PETSC_VIEWER_STDOUT_SELF` - standard output (default) 568 - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard 569 output where only the first processor opens 570 the file. All other processors send their 571 data to the first processor to print. 572 573 .seealso: [](ch_tao), `Tao`, `PetscViewerASCIIOpen()` 574 @*/ 575 PetscErrorCode TaoView(Tao tao, PetscViewer viewer) 576 { 577 PetscBool isascii, isstring; 578 TaoType type; 579 580 PetscFunctionBegin; 581 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 582 if (!viewer) PetscCall(PetscViewerASCIIGetStdout(((PetscObject)tao)->comm, &viewer)); 583 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 584 PetscCheckSameComm(tao, 1, viewer, 2); 585 586 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 587 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 588 if (isascii) { 589 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)tao, viewer)); 590 591 if (tao->ops->view) { 592 PetscCall(PetscViewerASCIIPushTab(viewer)); 593 PetscUseTypeMethod(tao, view, viewer); 594 PetscCall(PetscViewerASCIIPopTab(viewer)); 595 } 596 if (tao->linesearch) { 597 PetscCall(PetscViewerASCIIPushTab(viewer)); 598 PetscCall(TaoLineSearchView(tao->linesearch, viewer)); 599 PetscCall(PetscViewerASCIIPopTab(viewer)); 600 } 601 if (tao->ksp) { 602 PetscCall(PetscViewerASCIIPushTab(viewer)); 603 PetscCall(KSPView(tao->ksp, viewer)); 604 PetscCall(PetscViewerASCIIPrintf(viewer, "total KSP iterations: %" PetscInt_FMT "\n", tao->ksp_tot_its)); 605 PetscCall(PetscViewerASCIIPopTab(viewer)); 606 } 607 608 PetscCall(PetscViewerASCIIPushTab(viewer)); 609 610 if (tao->XL || tao->XU) PetscCall(PetscViewerASCIIPrintf(viewer, "Active Set subset type: %s\n", TaoSubSetTypes[tao->subset_type])); 611 612 PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: gatol=%g,", (double)tao->gatol)); 613 PetscCall(PetscViewerASCIIPrintf(viewer, " steptol=%g,", (double)tao->steptol)); 614 PetscCall(PetscViewerASCIIPrintf(viewer, " gttol=%g\n", (double)tao->gttol)); 615 PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Function/Gradient:=%g\n", (double)tao->residual)); 616 617 if (tao->constrained) { 618 PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances:")); 619 PetscCall(PetscViewerASCIIPrintf(viewer, " catol=%g,", (double)tao->catol)); 620 PetscCall(PetscViewerASCIIPrintf(viewer, " crtol=%g\n", (double)tao->crtol)); 621 PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Constraints:=%g\n", (double)tao->cnorm)); 622 } 623 624 if (tao->trust < tao->steptol) { 625 PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: steptol=%g\n", (double)tao->steptol)); 626 PetscCall(PetscViewerASCIIPrintf(viewer, "Final trust region radius:=%g\n", (double)tao->trust)); 627 } 628 629 if (tao->fmin > -1.e25) PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: function minimum=%g\n", (double)tao->fmin)); 630 PetscCall(PetscViewerASCIIPrintf(viewer, "Objective value=%g\n", (double)tao->fc)); 631 632 PetscCall(PetscViewerASCIIPrintf(viewer, "total number of iterations=%" PetscInt_FMT ", ", tao->niter)); 633 PetscCall(PetscViewerASCIIPrintf(viewer, " (max: %" PetscInt_FMT ")\n", tao->max_it)); 634 635 if (tao->nfuncs > 0) { 636 PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function evaluations=%" PetscInt_FMT ",", tao->nfuncs)); 637 PetscCall(PetscViewerASCIIPrintf(viewer, " max: %" PetscInt_FMT "\n", tao->max_funcs)); 638 } 639 if (tao->ngrads > 0) { 640 PetscCall(PetscViewerASCIIPrintf(viewer, "total number of gradient evaluations=%" PetscInt_FMT ",", tao->ngrads)); 641 PetscCall(PetscViewerASCIIPrintf(viewer, " max: %" PetscInt_FMT "\n", tao->max_funcs)); 642 } 643 if (tao->nfuncgrads > 0) { 644 PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function/gradient evaluations=%" PetscInt_FMT ",", tao->nfuncgrads)); 645 PetscCall(PetscViewerASCIIPrintf(viewer, " (max: %" PetscInt_FMT ")\n", tao->max_funcs)); 646 } 647 if (tao->nhess > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Hessian evaluations=%" PetscInt_FMT "\n", tao->nhess)); 648 if (tao->nconstraints > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of constraint function evaluations=%" PetscInt_FMT "\n", tao->nconstraints)); 649 if (tao->njac > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Jacobian evaluations=%" PetscInt_FMT "\n", tao->njac)); 650 651 if (tao->reason > 0) { 652 PetscCall(PetscViewerASCIIPrintf(viewer, "Solution converged: ")); 653 switch (tao->reason) { 654 case TAO_CONVERGED_GATOL: 655 PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)|| <= gatol\n")); 656 break; 657 case TAO_CONVERGED_GRTOL: 658 PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)||/|f(X)| <= grtol\n")); 659 break; 660 case TAO_CONVERGED_GTTOL: 661 PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)||/||g(X0)|| <= gttol\n")); 662 break; 663 case TAO_CONVERGED_STEPTOL: 664 PetscCall(PetscViewerASCIIPrintf(viewer, " Steptol -- step size small\n")); 665 break; 666 case TAO_CONVERGED_MINF: 667 PetscCall(PetscViewerASCIIPrintf(viewer, " Minf -- f < fmin\n")); 668 break; 669 case TAO_CONVERGED_USER: 670 PetscCall(PetscViewerASCIIPrintf(viewer, " User Terminated\n")); 671 break; 672 default: 673 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 674 break; 675 } 676 } else { 677 PetscCall(PetscViewerASCIIPrintf(viewer, "Solver terminated: %d", tao->reason)); 678 switch (tao->reason) { 679 case TAO_DIVERGED_MAXITS: 680 PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum Iterations\n")); 681 break; 682 case TAO_DIVERGED_NAN: 683 PetscCall(PetscViewerASCIIPrintf(viewer, " NAN or Inf encountered\n")); 684 break; 685 case TAO_DIVERGED_MAXFCN: 686 PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum Function Evaluations\n")); 687 break; 688 case TAO_DIVERGED_LS_FAILURE: 689 PetscCall(PetscViewerASCIIPrintf(viewer, " Line Search Failure\n")); 690 break; 691 case TAO_DIVERGED_TR_REDUCTION: 692 PetscCall(PetscViewerASCIIPrintf(viewer, " Trust Region too small\n")); 693 break; 694 case TAO_DIVERGED_USER: 695 PetscCall(PetscViewerASCIIPrintf(viewer, " User Terminated\n")); 696 break; 697 default: 698 PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 699 break; 700 } 701 } 702 PetscCall(PetscViewerASCIIPopTab(viewer)); 703 } else if (isstring) { 704 PetscCall(TaoGetType(tao, &type)); 705 PetscCall(PetscViewerStringSPrintf(viewer, " %-3.3s", type)); 706 } 707 PetscFunctionReturn(PETSC_SUCCESS); 708 } 709 710 /*@ 711 TaoSetRecycleHistory - Sets the boolean flag to enable/disable re-using 712 iterate information from the previous `TaoSolve()`. This feature is disabled by 713 default. 714 715 Logically Collective 716 717 Input Parameters: 718 + tao - the `Tao` context 719 - recycle - boolean flag 720 721 Options Database Key: 722 . -tao_recycle_history <true,false> - reuse the history 723 724 Level: intermediate 725 726 Notes: 727 For conjugate gradient methods (`TAOBNCG`), this re-uses the latest search direction 728 from the previous `TaoSolve()` call when computing the first search direction in a 729 new solution. By default, CG methods set the first search direction to the 730 negative gradient. 731 732 For quasi-Newton family of methods (`TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`), this re-uses 733 the accumulated quasi-Newton Hessian approximation from the previous `TaoSolve()` 734 call. By default, QN family of methods reset the initial Hessian approximation to 735 the identity matrix. 736 737 For any other algorithm, this setting has no effect. 738 739 .seealso: [](ch_tao), `Tao`, `TaoGetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL` 740 @*/ 741 PetscErrorCode TaoSetRecycleHistory(Tao tao, PetscBool recycle) 742 { 743 PetscFunctionBegin; 744 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 745 PetscValidLogicalCollectiveBool(tao, recycle, 2); 746 tao->recycle = recycle; 747 PetscFunctionReturn(PETSC_SUCCESS); 748 } 749 750 /*@ 751 TaoGetRecycleHistory - Retrieve the boolean flag for re-using iterate information 752 from the previous `TaoSolve()`. This feature is disabled by default. 753 754 Logically Collective 755 756 Input Parameter: 757 . tao - the `Tao` context 758 759 Output Parameter: 760 . recycle - boolean flag 761 762 Level: intermediate 763 764 .seealso: [](ch_tao), `Tao`, `TaoSetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL` 765 @*/ 766 PetscErrorCode TaoGetRecycleHistory(Tao tao, PetscBool *recycle) 767 { 768 PetscFunctionBegin; 769 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 770 PetscAssertPointer(recycle, 2); 771 *recycle = tao->recycle; 772 PetscFunctionReturn(PETSC_SUCCESS); 773 } 774 775 /*@ 776 TaoSetTolerances - Sets parameters used in `TaoSolve()` convergence tests 777 778 Logically Collective 779 780 Input Parameters: 781 + tao - the `Tao` context 782 . gatol - stop if norm of gradient is less than this 783 . grtol - stop if relative norm of gradient is less than this 784 - gttol - stop if norm of gradient is reduced by this factor 785 786 Options Database Keys: 787 + -tao_gatol <gatol> - Sets gatol 788 . -tao_grtol <grtol> - Sets grtol 789 - -tao_gttol <gttol> - Sets gttol 790 791 Stopping Criteria\: 792 .vb 793 ||g(X)|| <= gatol 794 ||g(X)|| / |f(X)| <= grtol 795 ||g(X)|| / ||g(X0)|| <= gttol 796 .ve 797 798 Level: beginner 799 800 Note: 801 Use `PETSC_DEFAULT` to leave one or more tolerances unchanged. 802 803 .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()` 804 @*/ 805 PetscErrorCode TaoSetTolerances(Tao tao, PetscReal gatol, PetscReal grtol, PetscReal gttol) 806 { 807 PetscFunctionBegin; 808 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 809 PetscValidLogicalCollectiveReal(tao, gatol, 2); 810 PetscValidLogicalCollectiveReal(tao, grtol, 3); 811 PetscValidLogicalCollectiveReal(tao, gttol, 4); 812 813 if (gatol != (PetscReal)PETSC_DEFAULT) { 814 if (gatol < 0) { 815 PetscCall(PetscInfo(tao, "Tried to set negative gatol -- ignored.\n")); 816 } else { 817 tao->gatol = PetscMax(0, gatol); 818 tao->gatol_changed = PETSC_TRUE; 819 } 820 } 821 822 if (grtol != (PetscReal)PETSC_DEFAULT) { 823 if (grtol < 0) { 824 PetscCall(PetscInfo(tao, "Tried to set negative grtol -- ignored.\n")); 825 } else { 826 tao->grtol = PetscMax(0, grtol); 827 tao->grtol_changed = PETSC_TRUE; 828 } 829 } 830 831 if (gttol != (PetscReal)PETSC_DEFAULT) { 832 if (gttol < 0) { 833 PetscCall(PetscInfo(tao, "Tried to set negative gttol -- ignored.\n")); 834 } else { 835 tao->gttol = PetscMax(0, gttol); 836 tao->gttol_changed = PETSC_TRUE; 837 } 838 } 839 PetscFunctionReturn(PETSC_SUCCESS); 840 } 841 842 /*@ 843 TaoSetConstraintTolerances - Sets constraint tolerance parameters used in `TaoSolve()` convergence tests 844 845 Logically Collective 846 847 Input Parameters: 848 + tao - the `Tao` context 849 . catol - absolute constraint tolerance, constraint norm must be less than `catol` for used for gatol convergence criteria 850 - crtol - relative constraint tolerance, constraint norm must be less than `crtol` for used for gatol, gttol convergence criteria 851 852 Options Database Keys: 853 + -tao_catol <catol> - Sets catol 854 - -tao_crtol <crtol> - Sets crtol 855 856 Level: intermediate 857 858 Notes: 859 Use `PETSC_DEFAULT` to leave any tolerance unchanged. 860 861 .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`, `TaoGetConstraintTolerances()`, `TaoSetTolerances()` 862 @*/ 863 PetscErrorCode TaoSetConstraintTolerances(Tao tao, PetscReal catol, PetscReal crtol) 864 { 865 PetscFunctionBegin; 866 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 867 PetscValidLogicalCollectiveReal(tao, catol, 2); 868 PetscValidLogicalCollectiveReal(tao, crtol, 3); 869 870 if (catol != (PetscReal)PETSC_DEFAULT) { 871 if (catol < 0) { 872 PetscCall(PetscInfo(tao, "Tried to set negative catol -- ignored.\n")); 873 } else { 874 tao->catol = PetscMax(0, catol); 875 tao->catol_changed = PETSC_TRUE; 876 } 877 } 878 879 if (crtol != (PetscReal)PETSC_DEFAULT) { 880 if (crtol < 0) { 881 PetscCall(PetscInfo(tao, "Tried to set negative crtol -- ignored.\n")); 882 } else { 883 tao->crtol = PetscMax(0, crtol); 884 tao->crtol_changed = PETSC_TRUE; 885 } 886 } 887 PetscFunctionReturn(PETSC_SUCCESS); 888 } 889 890 /*@ 891 TaoGetConstraintTolerances - Gets constraint tolerance parameters used in `TaoSolve()` convergence tests 892 893 Not Collective 894 895 Input Parameter: 896 . tao - the `Tao` context 897 898 Output Parameters: 899 + catol - absolute constraint tolerance, constraint norm must be less than `catol` for used for gatol convergence criteria 900 - crtol - relative constraint tolerance, constraint norm must be less than `crtol` for used for gatol, gttol convergence criteria 901 902 Level: intermediate 903 904 .seealso: [](ch_tao), `Tao`, `TaoConvergedReasons`,`TaoGetTolerances()`, `TaoSetTolerances()`, `TaoSetConstraintTolerances()` 905 @*/ 906 PetscErrorCode TaoGetConstraintTolerances(Tao tao, PetscReal *catol, PetscReal *crtol) 907 { 908 PetscFunctionBegin; 909 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 910 if (catol) *catol = tao->catol; 911 if (crtol) *crtol = tao->crtol; 912 PetscFunctionReturn(PETSC_SUCCESS); 913 } 914 915 /*@ 916 TaoSetFunctionLowerBound - Sets a bound on the solution objective value. 917 When an approximate solution with an objective value below this number 918 has been found, the solver will terminate. 919 920 Logically Collective 921 922 Input Parameters: 923 + tao - the Tao solver context 924 - fmin - the tolerance 925 926 Options Database Key: 927 . -tao_fmin <fmin> - sets the minimum function value 928 929 Level: intermediate 930 931 .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetTolerances()` 932 @*/ 933 PetscErrorCode TaoSetFunctionLowerBound(Tao tao, PetscReal fmin) 934 { 935 PetscFunctionBegin; 936 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 937 PetscValidLogicalCollectiveReal(tao, fmin, 2); 938 tao->fmin = fmin; 939 tao->fmin_changed = PETSC_TRUE; 940 PetscFunctionReturn(PETSC_SUCCESS); 941 } 942 943 /*@ 944 TaoGetFunctionLowerBound - Gets the bound on the solution objective value. 945 When an approximate solution with an objective value below this number 946 has been found, the solver will terminate. 947 948 Not Collective 949 950 Input Parameter: 951 . tao - the `Tao` solver context 952 953 Output Parameter: 954 . fmin - the minimum function value 955 956 Level: intermediate 957 958 .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetFunctionLowerBound()` 959 @*/ 960 PetscErrorCode TaoGetFunctionLowerBound(Tao tao, PetscReal *fmin) 961 { 962 PetscFunctionBegin; 963 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 964 PetscAssertPointer(fmin, 2); 965 *fmin = tao->fmin; 966 PetscFunctionReturn(PETSC_SUCCESS); 967 } 968 969 /*@ 970 TaoSetMaximumFunctionEvaluations - Sets a maximum number of function evaluations allowed for a `TaoSolve()`. 971 972 Logically Collective 973 974 Input Parameters: 975 + tao - the `Tao` solver context 976 - nfcn - the maximum number of function evaluations (>=0) 977 978 Options Database Key: 979 . -tao_max_funcs <nfcn> - sets the maximum number of function evaluations 980 981 Level: intermediate 982 983 .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoSetMaximumIterations()` 984 @*/ 985 PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao tao, PetscInt nfcn) 986 { 987 PetscFunctionBegin; 988 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 989 PetscValidLogicalCollectiveInt(tao, nfcn, 2); 990 if (nfcn >= 0) { 991 tao->max_funcs = PetscMax(0, nfcn); 992 } else { 993 tao->max_funcs = -1; 994 } 995 tao->max_funcs_changed = PETSC_TRUE; 996 PetscFunctionReturn(PETSC_SUCCESS); 997 } 998 999 /*@ 1000 TaoGetMaximumFunctionEvaluations - Gets a maximum number of function evaluations allowed for a `TaoSolve()` 1001 1002 Logically Collective 1003 1004 Input Parameter: 1005 . tao - the `Tao` solver context 1006 1007 Output Parameter: 1008 . nfcn - the maximum number of function evaluations 1009 1010 Level: intermediate 1011 1012 .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()` 1013 @*/ 1014 PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao tao, PetscInt *nfcn) 1015 { 1016 PetscFunctionBegin; 1017 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1018 PetscAssertPointer(nfcn, 2); 1019 *nfcn = tao->max_funcs; 1020 PetscFunctionReturn(PETSC_SUCCESS); 1021 } 1022 1023 /*@ 1024 TaoGetCurrentFunctionEvaluations - Get current number of function evaluations used by a `Tao` object 1025 1026 Not Collective 1027 1028 Input Parameter: 1029 . tao - the `Tao` solver context 1030 1031 Output Parameter: 1032 . nfuncs - the current number of function evaluations (maximum between gradient and function evaluations) 1033 1034 Level: intermediate 1035 1036 .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()` 1037 @*/ 1038 PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao tao, PetscInt *nfuncs) 1039 { 1040 PetscFunctionBegin; 1041 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1042 PetscAssertPointer(nfuncs, 2); 1043 *nfuncs = PetscMax(tao->nfuncs, tao->nfuncgrads); 1044 PetscFunctionReturn(PETSC_SUCCESS); 1045 } 1046 1047 /*@ 1048 TaoSetMaximumIterations - Sets a maximum number of iterates to be used in `TaoSolve()` 1049 1050 Logically Collective 1051 1052 Input Parameters: 1053 + tao - the `Tao` solver context 1054 - maxits - the maximum number of iterates (>=0) 1055 1056 Options Database Key: 1057 . -tao_max_it <its> - sets the maximum number of iterations 1058 1059 Level: intermediate 1060 1061 .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoSetMaximumFunctionEvaluations()` 1062 @*/ 1063 PetscErrorCode TaoSetMaximumIterations(Tao tao, PetscInt maxits) 1064 { 1065 PetscFunctionBegin; 1066 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1067 PetscValidLogicalCollectiveInt(tao, maxits, 2); 1068 tao->max_it = PetscMax(0, maxits); 1069 tao->max_it_changed = PETSC_TRUE; 1070 PetscFunctionReturn(PETSC_SUCCESS); 1071 } 1072 1073 /*@ 1074 TaoGetMaximumIterations - Gets a maximum number of iterates that will be used 1075 1076 Not Collective 1077 1078 Input Parameter: 1079 . tao - the `Tao` solver context 1080 1081 Output Parameter: 1082 . maxits - the maximum number of iterates 1083 1084 Level: intermediate 1085 1086 .seealso: [](ch_tao), `Tao`, `TaoSetMaximumIterations()`, `TaoGetMaximumFunctionEvaluations()` 1087 @*/ 1088 PetscErrorCode TaoGetMaximumIterations(Tao tao, PetscInt *maxits) 1089 { 1090 PetscFunctionBegin; 1091 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1092 PetscAssertPointer(maxits, 2); 1093 *maxits = tao->max_it; 1094 PetscFunctionReturn(PETSC_SUCCESS); 1095 } 1096 1097 /*@ 1098 TaoSetInitialTrustRegionRadius - Sets the initial trust region radius. 1099 1100 Logically Collective 1101 1102 Input Parameters: 1103 + tao - a `Tao` optimization solver 1104 - radius - the trust region radius 1105 1106 Options Database Key: 1107 . -tao_trust0 <t0> - sets initial trust region radius 1108 1109 Level: intermediate 1110 1111 .seealso: [](ch_tao), `Tao`, `TaoGetTrustRegionRadius()`, `TaoSetTrustRegionTolerance()`, `TAONTR` 1112 @*/ 1113 PetscErrorCode TaoSetInitialTrustRegionRadius(Tao tao, PetscReal radius) 1114 { 1115 PetscFunctionBegin; 1116 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1117 PetscValidLogicalCollectiveReal(tao, radius, 2); 1118 tao->trust0 = PetscMax(0.0, radius); 1119 tao->trust0_changed = PETSC_TRUE; 1120 PetscFunctionReturn(PETSC_SUCCESS); 1121 } 1122 1123 /*@ 1124 TaoGetInitialTrustRegionRadius - Gets the initial trust region radius. 1125 1126 Not Collective 1127 1128 Input Parameter: 1129 . tao - a `Tao` optimization solver 1130 1131 Output Parameter: 1132 . radius - the trust region radius 1133 1134 Level: intermediate 1135 1136 .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetCurrentTrustRegionRadius()`, `TAONTR` 1137 @*/ 1138 PetscErrorCode TaoGetInitialTrustRegionRadius(Tao tao, PetscReal *radius) 1139 { 1140 PetscFunctionBegin; 1141 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1142 PetscAssertPointer(radius, 2); 1143 *radius = tao->trust0; 1144 PetscFunctionReturn(PETSC_SUCCESS); 1145 } 1146 1147 /*@ 1148 TaoGetCurrentTrustRegionRadius - Gets the current trust region radius. 1149 1150 Not Collective 1151 1152 Input Parameter: 1153 . tao - a `Tao` optimization solver 1154 1155 Output Parameter: 1156 . radius - the trust region radius 1157 1158 Level: intermediate 1159 1160 .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetInitialTrustRegionRadius()`, `TAONTR` 1161 @*/ 1162 PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao tao, PetscReal *radius) 1163 { 1164 PetscFunctionBegin; 1165 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1166 PetscAssertPointer(radius, 2); 1167 *radius = tao->trust; 1168 PetscFunctionReturn(PETSC_SUCCESS); 1169 } 1170 1171 /*@ 1172 TaoGetTolerances - gets the current values of some tolerances used for the convergence testing of `TaoSolve()` 1173 1174 Not Collective 1175 1176 Input Parameter: 1177 . tao - the `Tao` context 1178 1179 Output Parameters: 1180 + gatol - stop if norm of gradient is less than this 1181 . grtol - stop if relative norm of gradient is less than this 1182 - gttol - stop if norm of gradient is reduced by a this factor 1183 1184 Level: intermediate 1185 1186 Note: 1187 `NULL` can be used as an argument if not all tolerances values are needed 1188 1189 .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()` 1190 @*/ 1191 PetscErrorCode TaoGetTolerances(Tao tao, PetscReal *gatol, PetscReal *grtol, PetscReal *gttol) 1192 { 1193 PetscFunctionBegin; 1194 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1195 if (gatol) *gatol = tao->gatol; 1196 if (grtol) *grtol = tao->grtol; 1197 if (gttol) *gttol = tao->gttol; 1198 PetscFunctionReturn(PETSC_SUCCESS); 1199 } 1200 1201 /*@ 1202 TaoGetKSP - Gets the linear solver used by the optimization solver. 1203 1204 Not Collective 1205 1206 Input Parameter: 1207 . tao - the `Tao` solver 1208 1209 Output Parameter: 1210 . ksp - the `KSP` linear solver used in the optimization solver 1211 1212 Level: intermediate 1213 1214 .seealso: [](ch_tao), `Tao`, `KSP` 1215 @*/ 1216 PetscErrorCode TaoGetKSP(Tao tao, KSP *ksp) 1217 { 1218 PetscFunctionBegin; 1219 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1220 PetscAssertPointer(ksp, 2); 1221 *ksp = tao->ksp; 1222 PetscFunctionReturn(PETSC_SUCCESS); 1223 } 1224 1225 /*@ 1226 TaoGetLinearSolveIterations - Gets the total number of linear iterations 1227 used by the `Tao` solver 1228 1229 Not Collective 1230 1231 Input Parameter: 1232 . tao - the `Tao` context 1233 1234 Output Parameter: 1235 . lits - number of linear iterations 1236 1237 Level: intermediate 1238 1239 Note: 1240 This counter is reset to zero for each successive call to `TaoSolve()` 1241 1242 .seealso: [](ch_tao), `Tao`, `TaoGetKSP()` 1243 @*/ 1244 PetscErrorCode TaoGetLinearSolveIterations(Tao tao, PetscInt *lits) 1245 { 1246 PetscFunctionBegin; 1247 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1248 PetscAssertPointer(lits, 2); 1249 *lits = tao->ksp_tot_its; 1250 PetscFunctionReturn(PETSC_SUCCESS); 1251 } 1252 1253 /*@ 1254 TaoGetLineSearch - Gets the line search used by the optimization solver. 1255 1256 Not Collective 1257 1258 Input Parameter: 1259 . tao - the `Tao` solver 1260 1261 Output Parameter: 1262 . ls - the line search used in the optimization solver 1263 1264 Level: intermediate 1265 1266 .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType` 1267 @*/ 1268 PetscErrorCode TaoGetLineSearch(Tao tao, TaoLineSearch *ls) 1269 { 1270 PetscFunctionBegin; 1271 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1272 PetscAssertPointer(ls, 2); 1273 *ls = tao->linesearch; 1274 PetscFunctionReturn(PETSC_SUCCESS); 1275 } 1276 1277 /*@ 1278 TaoAddLineSearchCounts - Adds the number of function evaluations spent 1279 in the line search to the running total. 1280 1281 Input Parameters: 1282 . tao - the `Tao` solver 1283 1284 Level: developer 1285 1286 .seealso: [](ch_tao), `Tao`, `TaoGetLineSearch()`, `TaoLineSearchApply()` 1287 @*/ 1288 PetscErrorCode TaoAddLineSearchCounts(Tao tao) 1289 { 1290 PetscBool flg; 1291 PetscInt nfeval, ngeval, nfgeval; 1292 1293 PetscFunctionBegin; 1294 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1295 if (tao->linesearch) { 1296 PetscCall(TaoLineSearchIsUsingTaoRoutines(tao->linesearch, &flg)); 1297 if (!flg) { 1298 PetscCall(TaoLineSearchGetNumberFunctionEvaluations(tao->linesearch, &nfeval, &ngeval, &nfgeval)); 1299 tao->nfuncs += nfeval; 1300 tao->ngrads += ngeval; 1301 tao->nfuncgrads += nfgeval; 1302 } 1303 } 1304 PetscFunctionReturn(PETSC_SUCCESS); 1305 } 1306 1307 /*@ 1308 TaoGetSolution - Returns the vector with the current solution from the `Tao` object 1309 1310 Not Collective 1311 1312 Input Parameter: 1313 . tao - the `Tao` context 1314 1315 Output Parameter: 1316 . X - the current solution 1317 1318 Level: intermediate 1319 1320 Note: 1321 The returned vector will be the same object that was passed into `TaoSetSolution()` 1322 1323 .seealso: [](ch_tao), `Tao`, `TaoSetSolution()`, `TaoSolve()` 1324 @*/ 1325 PetscErrorCode TaoGetSolution(Tao tao, Vec *X) 1326 { 1327 PetscFunctionBegin; 1328 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1329 PetscAssertPointer(X, 2); 1330 *X = tao->solution; 1331 PetscFunctionReturn(PETSC_SUCCESS); 1332 } 1333 1334 /*@ 1335 TaoResetStatistics - Initialize the statistics collected by the `Tao` object. 1336 These statistics include the iteration number, residual norms, and convergence status. 1337 This routine gets called before solving each optimization problem. 1338 1339 Collective 1340 1341 Input Parameter: 1342 . tao - the `Tao` context 1343 1344 Level: developer 1345 1346 .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()` 1347 @*/ 1348 PetscErrorCode TaoResetStatistics(Tao tao) 1349 { 1350 PetscFunctionBegin; 1351 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1352 tao->niter = 0; 1353 tao->nfuncs = 0; 1354 tao->nfuncgrads = 0; 1355 tao->ngrads = 0; 1356 tao->nhess = 0; 1357 tao->njac = 0; 1358 tao->nconstraints = 0; 1359 tao->ksp_its = 0; 1360 tao->ksp_tot_its = 0; 1361 tao->reason = TAO_CONTINUE_ITERATING; 1362 tao->residual = 0.0; 1363 tao->cnorm = 0.0; 1364 tao->step = 0.0; 1365 tao->lsflag = PETSC_FALSE; 1366 if (tao->hist_reset) tao->hist_len = 0; 1367 PetscFunctionReturn(PETSC_SUCCESS); 1368 } 1369 1370 /*@C 1371 TaoSetUpdate - Sets the general-purpose update function called 1372 at the beginning of every iteration of the optimization algorithm. Called after the new solution and the gradient 1373 is determined, but before the Hessian is computed (if applicable). 1374 1375 Logically Collective 1376 1377 Input Parameters: 1378 + tao - The `Tao` solver context 1379 - func - The function 1380 1381 Calling sequence of `func`: 1382 + tao - the optimizer context 1383 - ctx - The current step of the iteration 1384 1385 Level: advanced 1386 1387 .seealso: [](ch_tao), `Tao`, `TaoSolve()` 1388 @*/ 1389 PetscErrorCode TaoSetUpdate(Tao tao, PetscErrorCode (*func)(Tao, PetscInt, void *), void *ctx) 1390 { 1391 PetscFunctionBegin; 1392 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1393 tao->ops->update = func; 1394 tao->user_update = ctx; 1395 PetscFunctionReturn(PETSC_SUCCESS); 1396 } 1397 1398 /*@C 1399 TaoSetConvergenceTest - Sets the function that is to be used to test 1400 for convergence o fthe iterative minimization solution. The new convergence 1401 testing routine will replace Tao's default convergence test. 1402 1403 Logically Collective 1404 1405 Input Parameters: 1406 + tao - the `Tao` object 1407 . conv - the routine to test for convergence 1408 - ctx - [optional] context for private data for the convergence routine 1409 (may be `NULL`) 1410 1411 Calling sequence of `conv`: 1412 + tao - the `Tao` object 1413 - ctx - [optional] convergence context 1414 1415 Level: advanced 1416 1417 Note: 1418 The new convergence testing routine should call `TaoSetConvergedReason()`. 1419 1420 .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetConvergedReason()`, `TaoGetSolutionStatus()`, `TaoGetTolerances()`, `TaoSetMonitor` 1421 @*/ 1422 PetscErrorCode TaoSetConvergenceTest(Tao tao, PetscErrorCode (*conv)(Tao, void *), void *ctx) 1423 { 1424 PetscFunctionBegin; 1425 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1426 tao->ops->convergencetest = conv; 1427 tao->cnvP = ctx; 1428 PetscFunctionReturn(PETSC_SUCCESS); 1429 } 1430 1431 /*@C 1432 TaoSetMonitor - Sets an additional function that is to be used at every 1433 iteration of the solver to display the iteration's 1434 progress. 1435 1436 Logically Collective 1437 1438 Input Parameters: 1439 + tao - the `Tao` solver context 1440 . func - monitoring routine 1441 . ctx - [optional] user-defined context for private data for the monitor routine (may be `NULL`) 1442 - dest - [optional] function to destroy the context when the `Tao` is destroyed 1443 1444 Calling sequence of `func`: 1445 + tao - the `Tao` solver context 1446 - ctx - [optional] monitoring context 1447 1448 Calling sequence of `dest`: 1449 . ctx - monitoring context 1450 1451 Options Database Keys: 1452 + -tao_monitor - sets the default monitor `TaoMonitorDefault()` 1453 . -tao_smonitor - sets short monitor 1454 . -tao_cmonitor - same as smonitor plus constraint norm 1455 . -tao_view_solution - view solution at each iteration 1456 . -tao_view_gradient - view gradient at each iteration 1457 . -tao_view_ls_residual - view least-squares residual vector at each iteration 1458 - -tao_cancelmonitors - cancels all monitors that have been hardwired into a code by calls to TaoSetMonitor(), but does not cancel those set via the options database. 1459 1460 Level: intermediate 1461 1462 Notes: 1463 Several different monitoring routines may be set by calling 1464 `TaoSetMonitor()` multiple times; all will be called in the 1465 order in which they were set. 1466 1467 Fortran Notes: 1468 Only one monitor function may be set 1469 1470 .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoMonitorDefault()`, `TaoCancelMonitors()`, `TaoSetDestroyRoutine()`, `TaoView()` 1471 @*/ 1472 PetscErrorCode TaoSetMonitor(Tao tao, PetscErrorCode (*func)(Tao, void *), void *ctx, PetscErrorCode (*dest)(void **)) 1473 { 1474 PetscInt i; 1475 PetscBool identical; 1476 1477 PetscFunctionBegin; 1478 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1479 PetscCheck(tao->numbermonitors < MAXTAOMONITORS, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Cannot attach another monitor -- max=%d", MAXTAOMONITORS); 1480 1481 for (i = 0; i < tao->numbermonitors; i++) { 1482 PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))func, ctx, dest, (PetscErrorCode(*)(void))tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i], &identical)); 1483 if (identical) PetscFunctionReturn(PETSC_SUCCESS); 1484 } 1485 tao->monitor[tao->numbermonitors] = func; 1486 tao->monitorcontext[tao->numbermonitors] = (void *)ctx; 1487 tao->monitordestroy[tao->numbermonitors] = dest; 1488 ++tao->numbermonitors; 1489 PetscFunctionReturn(PETSC_SUCCESS); 1490 } 1491 1492 /*@ 1493 TaoCancelMonitors - Clears all the monitor functions for a `Tao` object. 1494 1495 Logically Collective 1496 1497 Input Parameter: 1498 . tao - the `Tao` solver context 1499 1500 Options Database Key: 1501 . -tao_cancelmonitors - cancels all monitors that have been hardwired 1502 into a code by calls to `TaoSetMonitor()`, but does not cancel those 1503 set via the options database 1504 1505 Level: advanced 1506 1507 Note: 1508 There is no way to clear one specific monitor from a `Tao` object. 1509 1510 .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoSetMonitor()` 1511 @*/ 1512 PetscErrorCode TaoCancelMonitors(Tao tao) 1513 { 1514 PetscInt i; 1515 1516 PetscFunctionBegin; 1517 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1518 for (i = 0; i < tao->numbermonitors; i++) { 1519 if (tao->monitordestroy[i]) PetscCall((*tao->monitordestroy[i])(&tao->monitorcontext[i])); 1520 } 1521 tao->numbermonitors = 0; 1522 PetscFunctionReturn(PETSC_SUCCESS); 1523 } 1524 1525 /*@ 1526 TaoMonitorDefault - Default routine for monitoring progress of `TaoSolve()` 1527 1528 Collective 1529 1530 Input Parameters: 1531 + tao - the `Tao` context 1532 - ctx - `PetscViewer` context or `NULL` 1533 1534 Options Database Key: 1535 . -tao_monitor - turn on default monitoring 1536 1537 Level: advanced 1538 1539 Note: 1540 This monitor prints the function value and gradient 1541 norm at each iteration. 1542 1543 .seealso: [](ch_tao), `Tao`, `TaoDefaultSMonitor()`, `TaoSetMonitor()` 1544 @*/ 1545 PetscErrorCode TaoMonitorDefault(Tao tao, void *ctx) 1546 { 1547 PetscInt its, tabs; 1548 PetscReal fct, gnorm; 1549 PetscViewer viewer = (PetscViewer)ctx; 1550 1551 PetscFunctionBegin; 1552 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1553 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1554 its = tao->niter; 1555 fct = tao->fc; 1556 gnorm = tao->residual; 1557 PetscCall(PetscViewerASCIIGetTab(viewer, &tabs)); 1558 PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel)); 1559 if (its == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) { 1560 PetscCall(PetscViewerASCIIPrintf(viewer, " Iteration information for %s solve.\n", ((PetscObject)tao)->prefix)); 1561 tao->header_printed = PETSC_TRUE; 1562 } 1563 PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " TAO,", its)); 1564 PetscCall(PetscViewerASCIIPrintf(viewer, " Function value: %g,", (double)fct)); 1565 if (gnorm >= PETSC_INFINITY) { 1566 PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: Inf \n")); 1567 } else { 1568 PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: %g \n", (double)gnorm)); 1569 } 1570 PetscCall(PetscViewerASCIISetTab(viewer, tabs)); 1571 PetscFunctionReturn(PETSC_SUCCESS); 1572 } 1573 1574 /*@ 1575 TaoDefaultGMonitor - Default routine for monitoring progress of `TaoSolve()` with extra detail on the globalization method. 1576 1577 Collective 1578 1579 Input Parameters: 1580 + tao - the `Tao` context 1581 - ctx - `PetscViewer` context or `NULL` 1582 1583 Options Database Key: 1584 . -tao_gmonitor - turn on monitoring with globalization information 1585 1586 Level: advanced 1587 1588 Note: 1589 This monitor prints the function value and gradient norm at each 1590 iteration, as well as the step size and trust radius. Note that the 1591 step size and trust radius may be the same for some algorithms. 1592 1593 .seealso: [](ch_tao), `Tao`, `TaoDefaultSMonitor()`, `TaoSetMonitor()` 1594 @*/ 1595 PetscErrorCode TaoDefaultGMonitor(Tao tao, void *ctx) 1596 { 1597 PetscInt its, tabs; 1598 PetscReal fct, gnorm, stp, tr; 1599 PetscViewer viewer = (PetscViewer)ctx; 1600 1601 PetscFunctionBegin; 1602 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1603 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1604 its = tao->niter; 1605 fct = tao->fc; 1606 gnorm = tao->residual; 1607 stp = tao->step; 1608 tr = tao->trust; 1609 PetscCall(PetscViewerASCIIGetTab(viewer, &tabs)); 1610 PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel)); 1611 if (its == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) { 1612 PetscCall(PetscViewerASCIIPrintf(viewer, " Iteration information for %s solve.\n", ((PetscObject)tao)->prefix)); 1613 tao->header_printed = PETSC_TRUE; 1614 } 1615 PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " TAO,", its)); 1616 PetscCall(PetscViewerASCIIPrintf(viewer, " Function value: %g,", (double)fct)); 1617 if (gnorm >= PETSC_INFINITY) { 1618 PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: Inf,")); 1619 } else { 1620 PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: %g,", (double)gnorm)); 1621 } 1622 PetscCall(PetscViewerASCIIPrintf(viewer, " Step: %g, Trust: %g\n", (double)stp, (double)tr)); 1623 PetscCall(PetscViewerASCIISetTab(viewer, tabs)); 1624 PetscFunctionReturn(PETSC_SUCCESS); 1625 } 1626 1627 /*@ 1628 TaoDefaultSMonitor - Default routine for monitoring progress of `TaoSolve()` 1629 1630 Collective 1631 1632 Input Parameters: 1633 + tao - the `Tao` context 1634 - ctx - `PetscViewer` context of type `PETSCVIEWERASCII` 1635 1636 Options Database Key: 1637 . -tao_smonitor - turn on default short monitoring 1638 1639 Level: advanced 1640 1641 Note: 1642 Same as `TaoMonitorDefault()` except 1643 it prints fewer digits of the residual as the residual gets smaller. 1644 This is because the later digits are meaningless and are often 1645 different on different machines; by using this routine different 1646 machines will usually generate the same output. 1647 1648 .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoSetMonitor()` 1649 @*/ 1650 PetscErrorCode TaoDefaultSMonitor(Tao tao, void *ctx) 1651 { 1652 PetscInt its, tabs; 1653 PetscReal fct, gnorm; 1654 PetscViewer viewer = (PetscViewer)ctx; 1655 1656 PetscFunctionBegin; 1657 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1658 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1659 its = tao->niter; 1660 fct = tao->fc; 1661 gnorm = tao->residual; 1662 PetscCall(PetscViewerASCIIGetTab(viewer, &tabs)); 1663 PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel)); 1664 PetscCall(PetscViewerASCIIPrintf(viewer, "iter = %3" PetscInt_FMT ",", its)); 1665 PetscCall(PetscViewerASCIIPrintf(viewer, " Function value %g,", (double)fct)); 1666 if (gnorm >= PETSC_INFINITY) { 1667 PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: Inf \n")); 1668 } else if (gnorm > 1.e-6) { 1669 PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: %g \n", (double)gnorm)); 1670 } else if (gnorm > 1.e-11) { 1671 PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: < 1.0e-6 \n")); 1672 } else { 1673 PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: < 1.0e-11 \n")); 1674 } 1675 PetscCall(PetscViewerASCIISetTab(viewer, tabs)); 1676 PetscFunctionReturn(PETSC_SUCCESS); 1677 } 1678 1679 /*@ 1680 TaoDefaultCMonitor - same as `TaoMonitorDefault()` except 1681 it prints the norm of the constraint function. 1682 1683 Collective 1684 1685 Input Parameters: 1686 + tao - the `Tao` context 1687 - ctx - `PetscViewer` context or `NULL` 1688 1689 Options Database Key: 1690 . -tao_cmonitor - monitor the constraints 1691 1692 Level: advanced 1693 1694 .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoSetMonitor()` 1695 @*/ 1696 PetscErrorCode TaoDefaultCMonitor(Tao tao, void *ctx) 1697 { 1698 PetscInt its, tabs; 1699 PetscReal fct, gnorm; 1700 PetscViewer viewer = (PetscViewer)ctx; 1701 1702 PetscFunctionBegin; 1703 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1704 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1705 its = tao->niter; 1706 fct = tao->fc; 1707 gnorm = tao->residual; 1708 PetscCall(PetscViewerASCIIGetTab(viewer, &tabs)); 1709 PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel)); 1710 PetscCall(PetscViewerASCIIPrintf(viewer, "iter = %" PetscInt_FMT ",", its)); 1711 PetscCall(PetscViewerASCIIPrintf(viewer, " Function value: %g,", (double)fct)); 1712 PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: %g ", (double)gnorm)); 1713 PetscCall(PetscViewerASCIIPrintf(viewer, " Constraint: %g \n", (double)tao->cnorm)); 1714 PetscCall(PetscViewerASCIISetTab(viewer, tabs)); 1715 PetscFunctionReturn(PETSC_SUCCESS); 1716 } 1717 1718 /*@C 1719 TaoSolutionMonitor - Views the solution at each iteration of `TaoSolve()` 1720 1721 Collective 1722 1723 Input Parameters: 1724 + tao - the `Tao` context 1725 - ctx - `PetscViewer` context or `NULL` 1726 1727 Options Database Key: 1728 . -tao_view_solution - view the solution 1729 1730 Level: advanced 1731 1732 .seealso: [](ch_tao), `Tao`, `TaoDefaultSMonitor()`, `TaoSetMonitor()` 1733 @*/ 1734 PetscErrorCode TaoSolutionMonitor(Tao tao, void *ctx) 1735 { 1736 PetscViewer viewer = (PetscViewer)ctx; 1737 1738 PetscFunctionBegin; 1739 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1740 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1741 PetscCall(VecView(tao->solution, viewer)); 1742 PetscFunctionReturn(PETSC_SUCCESS); 1743 } 1744 1745 /*@C 1746 TaoGradientMonitor - Views the gradient at each iteration of `TaoSolve()` 1747 1748 Collective 1749 1750 Input Parameters: 1751 + tao - the `Tao` context 1752 - ctx - `PetscViewer` context or `NULL` 1753 1754 Options Database Key: 1755 . -tao_view_gradient - view the gradient at each iteration 1756 1757 Level: advanced 1758 1759 .seealso: [](ch_tao), `Tao`, `TaoDefaultSMonitor()`, `TaoSetMonitor()` 1760 @*/ 1761 PetscErrorCode TaoGradientMonitor(Tao tao, void *ctx) 1762 { 1763 PetscViewer viewer = (PetscViewer)ctx; 1764 1765 PetscFunctionBegin; 1766 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1767 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1768 PetscCall(VecView(tao->gradient, viewer)); 1769 PetscFunctionReturn(PETSC_SUCCESS); 1770 } 1771 1772 /*@C 1773 TaoStepDirectionMonitor - Views the step-direction at each iteration of `TaoSolve()` 1774 1775 Collective 1776 1777 Input Parameters: 1778 + tao - the `Tao` context 1779 - ctx - `PetscViewer` context or `NULL` 1780 1781 Options Database Key: 1782 . -tao_view_stepdirection - view the step direction vector at each iteration 1783 1784 Level: advanced 1785 1786 .seealso: [](ch_tao), `Tao`, `TaoDefaultSMonitor()`, `TaoSetMonitor()` 1787 @*/ 1788 PetscErrorCode TaoStepDirectionMonitor(Tao tao, void *ctx) 1789 { 1790 PetscViewer viewer = (PetscViewer)ctx; 1791 1792 PetscFunctionBegin; 1793 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1794 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1795 PetscCall(VecView(tao->stepdirection, viewer)); 1796 PetscFunctionReturn(PETSC_SUCCESS); 1797 } 1798 1799 /*@C 1800 TaoDrawSolutionMonitor - Plots the solution at each iteration of `TaoSolve()` 1801 1802 Collective 1803 1804 Input Parameters: 1805 + tao - the `Tao` context 1806 - ctx - `TaoMonitorDraw` context 1807 1808 Options Database Key: 1809 . -tao_draw_solution - draw the solution at each iteration 1810 1811 Level: advanced 1812 1813 .seealso: [](ch_tao), `Tao`, `TaoSolutionMonitor()`, `TaoSetMonitor()`, `TaoDrawGradientMonitor`, `TaoMonitorDraw` 1814 @*/ 1815 PetscErrorCode TaoDrawSolutionMonitor(Tao tao, void *ctx) 1816 { 1817 TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx; 1818 1819 PetscFunctionBegin; 1820 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1821 if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(PETSC_SUCCESS); 1822 PetscCall(VecView(tao->solution, ictx->viewer)); 1823 PetscFunctionReturn(PETSC_SUCCESS); 1824 } 1825 1826 /*@C 1827 TaoDrawGradientMonitor - Plots the gradient at each iteration of `TaoSolve()` 1828 1829 Collective 1830 1831 Input Parameters: 1832 + tao - the `Tao` context 1833 - ctx - `PetscViewer` context 1834 1835 Options Database Key: 1836 . -tao_draw_gradient - draw the gradient at each iteration 1837 1838 Level: advanced 1839 1840 .seealso: [](ch_tao), `Tao`, `TaoGradientMonitor()`, `TaoSetMonitor()`, `TaoDrawSolutionMonitor` 1841 @*/ 1842 PetscErrorCode TaoDrawGradientMonitor(Tao tao, void *ctx) 1843 { 1844 TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx; 1845 1846 PetscFunctionBegin; 1847 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1848 if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(PETSC_SUCCESS); 1849 PetscCall(VecView(tao->gradient, ictx->viewer)); 1850 PetscFunctionReturn(PETSC_SUCCESS); 1851 } 1852 1853 /*@C 1854 TaoDrawStepMonitor - Plots the step direction at each iteration of `TaoSolve()` 1855 1856 Collective 1857 1858 Input Parameters: 1859 + tao - the `Tao` context 1860 - ctx - the `PetscViewer` context 1861 1862 Options Database Key: 1863 . -tao_draw_step - draw the step direction at each iteration 1864 1865 Level: advanced 1866 1867 .seealso: [](ch_tao), `Tao`, `TaoSetMonitor()`, `TaoDrawSolutionMonitor` 1868 @*/ 1869 PetscErrorCode TaoDrawStepMonitor(Tao tao, void *ctx) 1870 { 1871 PetscViewer viewer = (PetscViewer)ctx; 1872 1873 PetscFunctionBegin; 1874 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1875 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1876 PetscCall(VecView(tao->stepdirection, viewer)); 1877 PetscFunctionReturn(PETSC_SUCCESS); 1878 } 1879 1880 /*@C 1881 TaoResidualMonitor - Views the least-squares residual at each iteration of `TaoSolve()` 1882 1883 Collective 1884 1885 Input Parameters: 1886 + tao - the `Tao` context 1887 - ctx - the `PetscViewer` context or `NULL` 1888 1889 Options Database Key: 1890 . -tao_view_ls_residual - view the residual at each iteration 1891 1892 Level: advanced 1893 1894 .seealso: [](ch_tao), `Tao`, `TaoDefaultSMonitor()`, `TaoSetMonitor()` 1895 @*/ 1896 PetscErrorCode TaoResidualMonitor(Tao tao, void *ctx) 1897 { 1898 PetscViewer viewer = (PetscViewer)ctx; 1899 1900 PetscFunctionBegin; 1901 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1902 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1903 PetscCall(VecView(tao->ls_res, viewer)); 1904 PetscFunctionReturn(PETSC_SUCCESS); 1905 } 1906 1907 /*@ 1908 TaoDefaultConvergenceTest - Determines whether the solver should continue iterating 1909 or terminate. 1910 1911 Collective 1912 1913 Input Parameters: 1914 + tao - the `Tao` context 1915 - dummy - unused dummy context 1916 1917 Level: developer 1918 1919 Notes: 1920 This routine checks the residual in the optimality conditions, the 1921 relative residual in the optimity conditions, the number of function 1922 evaluations, and the function value to test convergence. Some 1923 solvers may use different convergence routines. 1924 1925 .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoGetConvergedReason()`, `TaoSetConvergedReason()` 1926 @*/ 1927 PetscErrorCode TaoDefaultConvergenceTest(Tao tao, void *dummy) 1928 { 1929 PetscInt niter = tao->niter, nfuncs = PetscMax(tao->nfuncs, tao->nfuncgrads); 1930 PetscInt max_funcs = tao->max_funcs; 1931 PetscReal gnorm = tao->residual, gnorm0 = tao->gnorm0; 1932 PetscReal f = tao->fc, steptol = tao->steptol, trradius = tao->step; 1933 PetscReal gatol = tao->gatol, grtol = tao->grtol, gttol = tao->gttol; 1934 PetscReal catol = tao->catol, crtol = tao->crtol; 1935 PetscReal fmin = tao->fmin, cnorm = tao->cnorm; 1936 TaoConvergedReason reason = tao->reason; 1937 1938 PetscFunctionBegin; 1939 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 1940 if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS); 1941 1942 if (PetscIsInfOrNanReal(f)) { 1943 PetscCall(PetscInfo(tao, "Failed to converged, function value is Inf or NaN\n")); 1944 reason = TAO_DIVERGED_NAN; 1945 } else if (f <= fmin && cnorm <= catol) { 1946 PetscCall(PetscInfo(tao, "Converged due to function value %g < minimum function value %g\n", (double)f, (double)fmin)); 1947 reason = TAO_CONVERGED_MINF; 1948 } else if (gnorm <= gatol && cnorm <= catol) { 1949 PetscCall(PetscInfo(tao, "Converged due to residual norm ||g(X)||=%g < %g\n", (double)gnorm, (double)gatol)); 1950 reason = TAO_CONVERGED_GATOL; 1951 } else if (f != 0 && PetscAbsReal(gnorm / f) <= grtol && cnorm <= crtol) { 1952 PetscCall(PetscInfo(tao, "Converged due to residual ||g(X)||/|f(X)| =%g < %g\n", (double)(gnorm / f), (double)grtol)); 1953 reason = TAO_CONVERGED_GRTOL; 1954 } else if (gnorm0 != 0 && ((gttol == 0 && gnorm == 0) || gnorm / gnorm0 < gttol) && cnorm <= crtol) { 1955 PetscCall(PetscInfo(tao, "Converged due to relative residual norm ||g(X)||/||g(X0)|| = %g < %g\n", (double)(gnorm / gnorm0), (double)gttol)); 1956 reason = TAO_CONVERGED_GTTOL; 1957 } else if (max_funcs >= 0 && nfuncs > max_funcs) { 1958 PetscCall(PetscInfo(tao, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", nfuncs, max_funcs)); 1959 reason = TAO_DIVERGED_MAXFCN; 1960 } else if (tao->lsflag != 0) { 1961 PetscCall(PetscInfo(tao, "Tao Line Search failure.\n")); 1962 reason = TAO_DIVERGED_LS_FAILURE; 1963 } else if (trradius < steptol && niter > 0) { 1964 PetscCall(PetscInfo(tao, "Trust region/step size too small: %g < %g\n", (double)trradius, (double)steptol)); 1965 reason = TAO_CONVERGED_STEPTOL; 1966 } else if (niter >= tao->max_it) { 1967 PetscCall(PetscInfo(tao, "Exceeded maximum number of iterations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", niter, tao->max_it)); 1968 reason = TAO_DIVERGED_MAXITS; 1969 } else { 1970 reason = TAO_CONTINUE_ITERATING; 1971 } 1972 tao->reason = reason; 1973 PetscFunctionReturn(PETSC_SUCCESS); 1974 } 1975 1976 /*@C 1977 TaoSetOptionsPrefix - Sets the prefix used for searching for all 1978 Tao options in the database. 1979 1980 Logically Collective 1981 1982 Input Parameters: 1983 + tao - the `Tao` context 1984 - p - the prefix string to prepend to all Tao option requests 1985 1986 Level: advanced 1987 1988 Notes: 1989 A hyphen (-) must NOT be given at the beginning of the prefix name. 1990 The first character of all runtime options is AUTOMATICALLY the hyphen. 1991 1992 For example, to distinguish between the runtime options for two 1993 different Tao solvers, one could call 1994 .vb 1995 TaoSetOptionsPrefix(tao1,"sys1_") 1996 TaoSetOptionsPrefix(tao2,"sys2_") 1997 .ve 1998 1999 This would enable use of different options for each system, such as 2000 .vb 2001 -sys1_tao_method blmvm -sys1_tao_grtol 1.e-3 2002 -sys2_tao_method lmvm -sys2_tao_grtol 1.e-4 2003 .ve 2004 2005 .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoAppendOptionsPrefix()`, `TaoGetOptionsPrefix()` 2006 @*/ 2007 PetscErrorCode TaoSetOptionsPrefix(Tao tao, const char p[]) 2008 { 2009 PetscFunctionBegin; 2010 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2011 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tao, p)); 2012 if (tao->linesearch) PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, p)); 2013 if (tao->ksp) PetscCall(KSPSetOptionsPrefix(tao->ksp, p)); 2014 PetscFunctionReturn(PETSC_SUCCESS); 2015 } 2016 2017 /*@C 2018 TaoAppendOptionsPrefix - Appends to the prefix used for searching for all Tao options in the database. 2019 2020 Logically Collective 2021 2022 Input Parameters: 2023 + tao - the `Tao` solver context 2024 - p - the prefix string to prepend to all `Tao` option requests 2025 2026 Level: advanced 2027 2028 Note: 2029 A hyphen (-) must NOT be given at the beginning of the prefix name. 2030 The first character of all runtime options is automatically the hyphen. 2031 2032 .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoSetOptionsPrefix()`, `TaoGetOptionsPrefix()` 2033 @*/ 2034 PetscErrorCode TaoAppendOptionsPrefix(Tao tao, const char p[]) 2035 { 2036 PetscFunctionBegin; 2037 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2038 PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao, p)); 2039 if (tao->linesearch) PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao->linesearch, p)); 2040 if (tao->ksp) PetscCall(KSPAppendOptionsPrefix(tao->ksp, p)); 2041 PetscFunctionReturn(PETSC_SUCCESS); 2042 } 2043 2044 /*@C 2045 TaoGetOptionsPrefix - Gets the prefix used for searching for all 2046 Tao options in the database 2047 2048 Not Collective 2049 2050 Input Parameter: 2051 . tao - the `Tao` context 2052 2053 Output Parameter: 2054 . p - pointer to the prefix string used is returned 2055 2056 Fortran Notes: 2057 Pass in a string 'prefix' of sufficient length to hold the prefix. 2058 2059 Level: advanced 2060 2061 .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoSetOptionsPrefix()`, `TaoAppendOptionsPrefix()` 2062 @*/ 2063 PetscErrorCode TaoGetOptionsPrefix(Tao tao, const char *p[]) 2064 { 2065 PetscFunctionBegin; 2066 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2067 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)tao, p)); 2068 PetscFunctionReturn(PETSC_SUCCESS); 2069 } 2070 2071 /*@C 2072 TaoSetType - Sets the `TaoType` for the minimization solver. 2073 2074 Collective 2075 2076 Input Parameters: 2077 + tao - the `Tao` solver context 2078 - type - a known method 2079 2080 Options Database Key: 2081 . -tao_type <type> - Sets the method; use -help for a list 2082 of available methods (for instance, "-tao_type lmvm" or "-tao_type tron") 2083 2084 Level: intermediate 2085 2086 .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoGetType()`, `TaoType` 2087 @*/ 2088 PetscErrorCode TaoSetType(Tao tao, TaoType type) 2089 { 2090 PetscErrorCode (*create_xxx)(Tao); 2091 PetscBool issame; 2092 2093 PetscFunctionBegin; 2094 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2095 2096 PetscCall(PetscObjectTypeCompare((PetscObject)tao, type, &issame)); 2097 if (issame) PetscFunctionReturn(PETSC_SUCCESS); 2098 2099 PetscCall(PetscFunctionListFind(TaoList, type, (void (**)(void)) & create_xxx)); 2100 PetscCheck(create_xxx, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Tao type %s", type); 2101 2102 /* Destroy the existing solver information */ 2103 PetscTryTypeMethod(tao, destroy); 2104 PetscCall(KSPDestroy(&tao->ksp)); 2105 PetscCall(TaoLineSearchDestroy(&tao->linesearch)); 2106 tao->ops->setup = NULL; 2107 tao->ops->solve = NULL; 2108 tao->ops->view = NULL; 2109 tao->ops->setfromoptions = NULL; 2110 tao->ops->destroy = NULL; 2111 2112 tao->setupcalled = PETSC_FALSE; 2113 2114 PetscCall((*create_xxx)(tao)); 2115 PetscCall(PetscObjectChangeTypeName((PetscObject)tao, type)); 2116 PetscFunctionReturn(PETSC_SUCCESS); 2117 } 2118 2119 /*@C 2120 TaoRegister - Adds a method to the Tao package for minimization. 2121 2122 Not Collective 2123 2124 Input Parameters: 2125 + sname - name of a new user-defined solver 2126 - func - routine to Create method context 2127 2128 Example Usage: 2129 .vb 2130 TaoRegister("my_solver", MySolverCreate); 2131 .ve 2132 2133 Then, your solver can be chosen with the procedural interface via 2134 $ TaoSetType(tao, "my_solver") 2135 or at runtime via the option 2136 $ -tao_type my_solver 2137 2138 Level: advanced 2139 2140 Note: 2141 `TaoRegister()` may be called multiple times to add several user-defined solvers. 2142 2143 .seealso: [](ch_tao), `Tao`, `TaoSetType()`, `TaoRegisterAll()`, `TaoRegisterDestroy()` 2144 @*/ 2145 PetscErrorCode TaoRegister(const char sname[], PetscErrorCode (*func)(Tao)) 2146 { 2147 PetscFunctionBegin; 2148 PetscCall(TaoInitializePackage()); 2149 PetscCall(PetscFunctionListAdd(&TaoList, sname, (void (*)(void))func)); 2150 PetscFunctionReturn(PETSC_SUCCESS); 2151 } 2152 2153 /*@C 2154 TaoRegisterDestroy - Frees the list of minimization solvers that were 2155 registered by `TaoRegister()`. 2156 2157 Not Collective 2158 2159 Level: advanced 2160 2161 .seealso: [](ch_tao), `Tao`, `TaoRegisterAll()`, `TaoRegister()` 2162 @*/ 2163 PetscErrorCode TaoRegisterDestroy(void) 2164 { 2165 PetscFunctionBegin; 2166 PetscCall(PetscFunctionListDestroy(&TaoList)); 2167 TaoRegisterAllCalled = PETSC_FALSE; 2168 PetscFunctionReturn(PETSC_SUCCESS); 2169 } 2170 2171 /*@ 2172 TaoGetIterationNumber - Gets the number of `TaoSolve()` iterations completed 2173 at this time. 2174 2175 Not Collective 2176 2177 Input Parameter: 2178 . tao - the `Tao` context 2179 2180 Output Parameter: 2181 . iter - iteration number 2182 2183 Notes: 2184 For example, during the computation of iteration 2 this would return 1. 2185 2186 Level: intermediate 2187 2188 .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`, `TaoGetResidualNorm()`, `TaoGetObjective()` 2189 @*/ 2190 PetscErrorCode TaoGetIterationNumber(Tao tao, PetscInt *iter) 2191 { 2192 PetscFunctionBegin; 2193 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2194 PetscAssertPointer(iter, 2); 2195 *iter = tao->niter; 2196 PetscFunctionReturn(PETSC_SUCCESS); 2197 } 2198 2199 /*@ 2200 TaoGetResidualNorm - Gets the current value of the norm of the residual (gradient) 2201 at this time. 2202 2203 Not Collective 2204 2205 Input Parameter: 2206 . tao - the `Tao` context 2207 2208 Output Parameter: 2209 . value - the current value 2210 2211 Level: intermediate 2212 2213 Developer Notes: 2214 This is the 2-norm of the residual, we cannot use `TaoGetGradientNorm()` because that has 2215 a different meaning. For some reason `Tao` sometimes calls the gradient the residual. 2216 2217 .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`, `TaoGetIterationNumber()`, `TaoGetObjective()` 2218 @*/ 2219 PetscErrorCode TaoGetResidualNorm(Tao tao, PetscReal *value) 2220 { 2221 PetscFunctionBegin; 2222 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2223 PetscAssertPointer(value, 2); 2224 *value = tao->residual; 2225 PetscFunctionReturn(PETSC_SUCCESS); 2226 } 2227 2228 /*@ 2229 TaoSetIterationNumber - Sets the current iteration number. 2230 2231 Logically Collective 2232 2233 Input Parameters: 2234 + tao - the `Tao` context 2235 - iter - iteration number 2236 2237 Level: developer 2238 2239 .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()` 2240 @*/ 2241 PetscErrorCode TaoSetIterationNumber(Tao tao, PetscInt iter) 2242 { 2243 PetscFunctionBegin; 2244 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2245 PetscValidLogicalCollectiveInt(tao, iter, 2); 2246 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao)); 2247 tao->niter = iter; 2248 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao)); 2249 PetscFunctionReturn(PETSC_SUCCESS); 2250 } 2251 2252 /*@ 2253 TaoGetTotalIterationNumber - Gets the total number of `TaoSolve()` iterations 2254 completed. This number keeps accumulating if multiple solves 2255 are called with the `Tao` object. 2256 2257 Not Collective 2258 2259 Input Parameter: 2260 . tao - the `Tao` context 2261 2262 Output Parameter: 2263 . iter - number of iterations 2264 2265 Level: intermediate 2266 2267 Note: 2268 The total iteration count is updated after each solve, if there is a current 2269 `TaoSolve()` in progress then those iterations are not included in the count 2270 2271 .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()` 2272 @*/ 2273 PetscErrorCode TaoGetTotalIterationNumber(Tao tao, PetscInt *iter) 2274 { 2275 PetscFunctionBegin; 2276 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2277 PetscAssertPointer(iter, 2); 2278 *iter = tao->ntotalits; 2279 PetscFunctionReturn(PETSC_SUCCESS); 2280 } 2281 2282 /*@ 2283 TaoSetTotalIterationNumber - Sets the current total iteration number. 2284 2285 Logically Collective 2286 2287 Input Parameters: 2288 + tao - the `Tao` context 2289 - iter - the iteration number 2290 2291 Level: developer 2292 2293 .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()` 2294 @*/ 2295 PetscErrorCode TaoSetTotalIterationNumber(Tao tao, PetscInt iter) 2296 { 2297 PetscFunctionBegin; 2298 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2299 PetscValidLogicalCollectiveInt(tao, iter, 2); 2300 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao)); 2301 tao->ntotalits = iter; 2302 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao)); 2303 PetscFunctionReturn(PETSC_SUCCESS); 2304 } 2305 2306 /*@ 2307 TaoSetConvergedReason - Sets the termination flag on a `Tao` object 2308 2309 Logically Collective 2310 2311 Input Parameters: 2312 + tao - the `Tao` context 2313 - reason - the `TaoConvergedReason` 2314 2315 Level: intermediate 2316 2317 .seealso: [](ch_tao), `Tao`, `TaoConvergedReason` 2318 @*/ 2319 PetscErrorCode TaoSetConvergedReason(Tao tao, TaoConvergedReason reason) 2320 { 2321 PetscFunctionBegin; 2322 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2323 PetscValidLogicalCollectiveEnum(tao, reason, 2); 2324 tao->reason = reason; 2325 PetscFunctionReturn(PETSC_SUCCESS); 2326 } 2327 2328 /*@ 2329 TaoGetConvergedReason - Gets the reason the `TaoSolve()` was stopped. 2330 2331 Not Collective 2332 2333 Input Parameter: 2334 . tao - the `Tao` solver context 2335 2336 Output Parameter: 2337 . reason - value of `TaoConvergedReason` 2338 2339 Level: intermediate 2340 2341 .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetConvergenceTest()`, `TaoSetTolerances()` 2342 @*/ 2343 PetscErrorCode TaoGetConvergedReason(Tao tao, TaoConvergedReason *reason) 2344 { 2345 PetscFunctionBegin; 2346 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2347 PetscAssertPointer(reason, 2); 2348 *reason = tao->reason; 2349 PetscFunctionReturn(PETSC_SUCCESS); 2350 } 2351 2352 /*@ 2353 TaoGetSolutionStatus - Get the current iterate, objective value, 2354 residual, infeasibility, and termination from a `Tao` object 2355 2356 Not Collective 2357 2358 Input Parameter: 2359 . tao - the `Tao` context 2360 2361 Output Parameters: 2362 + its - the current iterate number (>=0) 2363 . f - the current function value 2364 . gnorm - the square of the gradient norm, duality gap, or other measure indicating distance from optimality. 2365 . cnorm - the infeasibility of the current solution with regard to the constraints. 2366 . xdiff - the step length or trust region radius of the most recent iterate. 2367 - reason - The termination reason, which can equal `TAO_CONTINUE_ITERATING` 2368 2369 Level: intermediate 2370 2371 Notes: 2372 Tao returns the values set by the solvers in the routine `TaoMonitor()`. 2373 2374 If any of the output arguments are set to `NULL`, no corresponding value will be returned. 2375 2376 .seealso: [](ch_tao), `TaoMonitor()`, `TaoGetConvergedReason()` 2377 @*/ 2378 PetscErrorCode TaoGetSolutionStatus(Tao tao, PetscInt *its, PetscReal *f, PetscReal *gnorm, PetscReal *cnorm, PetscReal *xdiff, TaoConvergedReason *reason) 2379 { 2380 PetscFunctionBegin; 2381 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2382 if (its) *its = tao->niter; 2383 if (f) *f = tao->fc; 2384 if (gnorm) *gnorm = tao->residual; 2385 if (cnorm) *cnorm = tao->cnorm; 2386 if (reason) *reason = tao->reason; 2387 if (xdiff) *xdiff = tao->step; 2388 PetscFunctionReturn(PETSC_SUCCESS); 2389 } 2390 2391 /*@C 2392 TaoGetType - Gets the current `TaoType` being used in the `Tao` object 2393 2394 Not Collective 2395 2396 Input Parameter: 2397 . tao - the `Tao` solver context 2398 2399 Output Parameter: 2400 . type - the `TaoType` 2401 2402 Level: intermediate 2403 2404 .seealso: [](ch_tao), `Tao`, `TaoType`, `TaoSetType()` 2405 @*/ 2406 PetscErrorCode TaoGetType(Tao tao, TaoType *type) 2407 { 2408 PetscFunctionBegin; 2409 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2410 PetscAssertPointer(type, 2); 2411 *type = ((PetscObject)tao)->type_name; 2412 PetscFunctionReturn(PETSC_SUCCESS); 2413 } 2414 2415 /*@C 2416 TaoMonitor - Monitor the solver and the current solution. This 2417 routine will record the iteration number and residual statistics, 2418 and call any monitors specified by the user. 2419 2420 Input Parameters: 2421 + tao - the `Tao` context 2422 . its - the current iterate number (>=0) 2423 . f - the current objective function value 2424 . res - the gradient norm, square root of the duality gap, or other measure indicating distance from optimality. This measure will be recorded and 2425 used for some termination tests. 2426 . cnorm - the infeasibility of the current solution with regard to the constraints. 2427 - steplength - multiple of the step direction added to the previous iterate. 2428 2429 Options Database Key: 2430 . -tao_monitor - Use the default monitor, which prints statistics to standard output 2431 2432 Level: developer 2433 2434 .seealso: [](ch_tao), `Tao`, `TaoGetConvergedReason()`, `TaoMonitorDefault()`, `TaoSetMonitor()` 2435 @*/ 2436 PetscErrorCode TaoMonitor(Tao tao, PetscInt its, PetscReal f, PetscReal res, PetscReal cnorm, PetscReal steplength) 2437 { 2438 PetscInt i; 2439 2440 PetscFunctionBegin; 2441 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2442 tao->fc = f; 2443 tao->residual = res; 2444 tao->cnorm = cnorm; 2445 tao->step = steplength; 2446 if (!its) { 2447 tao->cnorm0 = cnorm; 2448 tao->gnorm0 = res; 2449 } 2450 PetscCall(VecLockReadPush(tao->solution)); 2451 for (i = 0; i < tao->numbermonitors; i++) PetscCall((*tao->monitor[i])(tao, tao->monitorcontext[i])); 2452 PetscCall(VecLockReadPop(tao->solution)); 2453 PetscFunctionReturn(PETSC_SUCCESS); 2454 } 2455 2456 /*@ 2457 TaoSetConvergenceHistory - Sets the array used to hold the convergence history. 2458 2459 Logically Collective 2460 2461 Input Parameters: 2462 + tao - the `Tao` solver context 2463 . obj - array to hold objective value history 2464 . resid - array to hold residual history 2465 . cnorm - array to hold constraint violation history 2466 . lits - integer array holds the number of linear iterations for each Tao iteration 2467 . na - size of `obj`, `resid`, and `cnorm` 2468 - reset - `PETSC_TRUE` indicates each new minimization resets the history counter to zero, 2469 else it continues storing new values for new minimizations after the old ones 2470 2471 Level: intermediate 2472 2473 Notes: 2474 If set, `Tao` will fill the given arrays with the indicated 2475 information at each iteration. If 'obj','resid','cnorm','lits' are 2476 *all* `NULL` then space (using size `na`, or 1000 if na is `PETSC_DECIDE` or 2477 `PETSC_DEFAULT`) is allocated for the history. 2478 If not all are `NULL`, then only the non-`NULL` information categories 2479 will be stored, the others will be ignored. 2480 2481 Any convergence information after iteration number 'na' will not be stored. 2482 2483 This routine is useful, e.g., when running a code for purposes 2484 of accurate performance monitoring, when no I/O should be done 2485 during the section of code that is being timed. 2486 2487 .seealso: [](ch_tao), `TaoGetConvergenceHistory()` 2488 @*/ 2489 PetscErrorCode TaoSetConvergenceHistory(Tao tao, PetscReal obj[], PetscReal resid[], PetscReal cnorm[], PetscInt lits[], PetscInt na, PetscBool reset) 2490 { 2491 PetscFunctionBegin; 2492 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2493 if (obj) PetscAssertPointer(obj, 2); 2494 if (resid) PetscAssertPointer(resid, 3); 2495 if (cnorm) PetscAssertPointer(cnorm, 4); 2496 if (lits) PetscAssertPointer(lits, 5); 2497 2498 if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000; 2499 if (!obj && !resid && !cnorm && !lits) { 2500 PetscCall(PetscCalloc4(na, &obj, na, &resid, na, &cnorm, na, &lits)); 2501 tao->hist_malloc = PETSC_TRUE; 2502 } 2503 2504 tao->hist_obj = obj; 2505 tao->hist_resid = resid; 2506 tao->hist_cnorm = cnorm; 2507 tao->hist_lits = lits; 2508 tao->hist_max = na; 2509 tao->hist_reset = reset; 2510 tao->hist_len = 0; 2511 PetscFunctionReturn(PETSC_SUCCESS); 2512 } 2513 2514 /*@C 2515 TaoGetConvergenceHistory - Gets the arrays used that hold the convergence history. 2516 2517 Collective 2518 2519 Input Parameter: 2520 . tao - the `Tao` context 2521 2522 Output Parameters: 2523 + obj - array used to hold objective value history 2524 . resid - array used to hold residual history 2525 . cnorm - array used to hold constraint violation history 2526 . lits - integer array used to hold linear solver iteration count 2527 - nhist - size of `obj`, `resid`, `cnorm`, and `lits` 2528 2529 Level: advanced 2530 2531 Notes: 2532 This routine must be preceded by calls to `TaoSetConvergenceHistory()` 2533 and `TaoSolve()`, otherwise it returns useless information. 2534 2535 This routine is useful, e.g., when running a code for purposes 2536 of accurate performance monitoring, when no I/O should be done 2537 during the section of code that is being timed. 2538 2539 Fortran Notes: 2540 The calling sequence is 2541 .vb 2542 call TaoGetConvergenceHistory(Tao tao, PetscInt nhist, PetscErrorCode ierr) 2543 .ve 2544 2545 .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetConvergenceHistory()` 2546 @*/ 2547 PetscErrorCode TaoGetConvergenceHistory(Tao tao, PetscReal **obj, PetscReal **resid, PetscReal **cnorm, PetscInt **lits, PetscInt *nhist) 2548 { 2549 PetscFunctionBegin; 2550 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2551 if (obj) *obj = tao->hist_obj; 2552 if (cnorm) *cnorm = tao->hist_cnorm; 2553 if (resid) *resid = tao->hist_resid; 2554 if (lits) *lits = tao->hist_lits; 2555 if (nhist) *nhist = tao->hist_len; 2556 PetscFunctionReturn(PETSC_SUCCESS); 2557 } 2558 2559 /*@ 2560 TaoSetApplicationContext - Sets the optional user-defined context for a `Tao` solver. 2561 2562 Logically Collective 2563 2564 Input Parameters: 2565 + tao - the `Tao` context 2566 - usrP - optional user context 2567 2568 Level: intermediate 2569 2570 .seealso: [](ch_tao), `Tao`, `TaoGetApplicationContext()` 2571 @*/ 2572 PetscErrorCode TaoSetApplicationContext(Tao tao, void *usrP) 2573 { 2574 PetscFunctionBegin; 2575 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2576 tao->user = usrP; 2577 PetscFunctionReturn(PETSC_SUCCESS); 2578 } 2579 2580 /*@ 2581 TaoGetApplicationContext - Gets the user-defined context for a `Tao` solver 2582 2583 Not Collective 2584 2585 Input Parameter: 2586 . tao - the `Tao` context 2587 2588 Output Parameter: 2589 . usrP - user context 2590 2591 Level: intermediate 2592 2593 .seealso: [](ch_tao), `Tao`, `TaoSetApplicationContext()` 2594 @*/ 2595 PetscErrorCode TaoGetApplicationContext(Tao tao, void *usrP) 2596 { 2597 PetscFunctionBegin; 2598 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2599 PetscAssertPointer(usrP, 2); 2600 *(void **)usrP = tao->user; 2601 PetscFunctionReturn(PETSC_SUCCESS); 2602 } 2603 2604 /*@ 2605 TaoSetGradientNorm - Sets the matrix used to define the norm that measures the size of the gradient. 2606 2607 Collective 2608 2609 Input Parameters: 2610 + tao - the `Tao` context 2611 - M - matrix that defines the norm 2612 2613 Level: beginner 2614 2615 .seealso: [](ch_tao), `Tao`, `TaoGetGradientNorm()`, `TaoGradientNorm()` 2616 @*/ 2617 PetscErrorCode TaoSetGradientNorm(Tao tao, Mat M) 2618 { 2619 PetscFunctionBegin; 2620 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2621 PetscValidHeaderSpecific(M, MAT_CLASSID, 2); 2622 PetscCall(PetscObjectReference((PetscObject)M)); 2623 PetscCall(MatDestroy(&tao->gradient_norm)); 2624 PetscCall(VecDestroy(&tao->gradient_norm_tmp)); 2625 tao->gradient_norm = M; 2626 PetscCall(MatCreateVecs(M, NULL, &tao->gradient_norm_tmp)); 2627 PetscFunctionReturn(PETSC_SUCCESS); 2628 } 2629 2630 /*@ 2631 TaoGetGradientNorm - Returns the matrix used to define the norm used for measuring the size of the gradient. 2632 2633 Not Collective 2634 2635 Input Parameter: 2636 . tao - the `Tao` context 2637 2638 Output Parameter: 2639 . M - gradient norm 2640 2641 Level: beginner 2642 2643 .seealso: [](ch_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGradientNorm()` 2644 @*/ 2645 PetscErrorCode TaoGetGradientNorm(Tao tao, Mat *M) 2646 { 2647 PetscFunctionBegin; 2648 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2649 PetscAssertPointer(M, 2); 2650 *M = tao->gradient_norm; 2651 PetscFunctionReturn(PETSC_SUCCESS); 2652 } 2653 2654 /*@C 2655 TaoGradientNorm - Compute the norm using the `NormType`, the user has selected 2656 2657 Collective 2658 2659 Input Parameters: 2660 + tao - the `Tao` context 2661 . gradient - the gradient to be computed 2662 - type - the norm type 2663 2664 Output Parameter: 2665 . gnorm - the gradient norm 2666 2667 Level: advanced 2668 2669 .seealso: [](ch_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGetGradientNorm()` 2670 @*/ 2671 PetscErrorCode TaoGradientNorm(Tao tao, Vec gradient, NormType type, PetscReal *gnorm) 2672 { 2673 PetscFunctionBegin; 2674 PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 2675 PetscValidHeaderSpecific(gradient, VEC_CLASSID, 2); 2676 PetscValidLogicalCollectiveEnum(tao, type, 3); 2677 PetscAssertPointer(gnorm, 4); 2678 if (tao->gradient_norm) { 2679 PetscScalar gnorms; 2680 2681 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."); 2682 PetscCall(MatMult(tao->gradient_norm, gradient, tao->gradient_norm_tmp)); 2683 PetscCall(VecDot(gradient, tao->gradient_norm_tmp, &gnorms)); 2684 *gnorm = PetscRealPart(PetscSqrtScalar(gnorms)); 2685 } else { 2686 PetscCall(VecNorm(gradient, type, gnorm)); 2687 } 2688 PetscFunctionReturn(PETSC_SUCCESS); 2689 } 2690 2691 /*@C 2692 TaoMonitorDrawCtxCreate - Creates the monitor context for `TaoMonitorDrawSolution()` 2693 2694 Collective 2695 2696 Input Parameters: 2697 + comm - the communicator to share the context 2698 . host - the name of the X Windows host that will display the monitor 2699 . label - the label to put at the top of the display window 2700 . x - the horizontal coordinate of the lower left corner of the window to open 2701 . y - the vertical coordinate of the lower left corner of the window to open 2702 . m - the width of the window 2703 . n - the height of the window 2704 - howoften - how many `Tao` iterations between displaying the monitor information 2705 2706 Output Parameter: 2707 . ctx - the monitor context 2708 2709 Options Database Key: 2710 . -tao_draw_solution_initial - show initial guess as well as current solution 2711 2712 Level: intermediate 2713 2714 .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorDrawCtx()` 2715 @*/ 2716 PetscErrorCode TaoMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TaoMonitorDrawCtx *ctx) 2717 { 2718 PetscFunctionBegin; 2719 PetscCall(PetscNew(ctx)); 2720 PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer)); 2721 PetscCall(PetscViewerSetFromOptions((*ctx)->viewer)); 2722 (*ctx)->howoften = howoften; 2723 PetscFunctionReturn(PETSC_SUCCESS); 2724 } 2725 2726 /*@C 2727 TaoMonitorDrawCtxDestroy - Destroys the monitor context for `TaoMonitorDrawSolution()` 2728 2729 Collective 2730 2731 Input Parameter: 2732 . ictx - the monitor context 2733 2734 Level: intermediate 2735 2736 .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorDrawSolution()` 2737 @*/ 2738 PetscErrorCode TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx *ictx) 2739 { 2740 PetscFunctionBegin; 2741 PetscCall(PetscViewerDestroy(&(*ictx)->viewer)); 2742 PetscCall(PetscFree(*ictx)); 2743 PetscFunctionReturn(PETSC_SUCCESS); 2744 } 2745