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