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