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