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