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