1 #include <petsc/private/linesearchimpl.h> /*I "petscsnes.h" I*/ 2 3 PetscBool SNESLineSearchRegisterAllCalled = PETSC_FALSE; 4 PetscFunctionList SNESLineSearchList = NULL; 5 6 PetscClassId SNESLINESEARCH_CLASSID; 7 PetscLogEvent SNESLINESEARCH_Apply; 8 9 /*@ 10 SNESLineSearchMonitorCancel - Clears all the monitor functions for a `SNESLineSearch` object. 11 12 Logically Collective 13 14 Input Parameter: 15 . ls - the `SNESLineSearch` context 16 17 Options Database Key: 18 . -snes_linesearch_monitor_cancel - cancels all monitors that have been hardwired 19 into a code by calls to `SNESLineSearchMonitorSet()`, but does not cancel those 20 set via the options database 21 22 Level: advanced 23 24 Notes: 25 There is no way to clear one specific monitor from a `SNESLineSearch` object. 26 27 This does not clear the monitor set with `SNESLineSearchSetDefaultMonitor()` use `SNESLineSearchSetDefaultMonitor`(ls,NULL) to cancel 28 that one. 29 30 .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorSet()` 31 @*/ 32 PetscErrorCode SNESLineSearchMonitorCancel(SNESLineSearch ls) 33 { 34 PetscInt i; 35 36 PetscFunctionBegin; 37 PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1); 38 for (i = 0; i < ls->numbermonitors; i++) { 39 if (ls->monitordestroy[i]) PetscCall((*ls->monitordestroy[i])(&ls->monitorcontext[i])); 40 } 41 ls->numbermonitors = 0; 42 PetscFunctionReturn(PETSC_SUCCESS); 43 } 44 45 /*@ 46 SNESLineSearchMonitor - runs the user provided monitor routines, if they exist 47 48 Collective 49 50 Input Parameter: 51 . ls - the linesearch object 52 53 Level: developer 54 55 Note: 56 This routine is called by the `SNES` implementations. 57 It does not typically need to be called by the user. 58 59 .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorSet()` 60 @*/ 61 PetscErrorCode SNESLineSearchMonitor(SNESLineSearch ls) 62 { 63 PetscInt i, n = ls->numbermonitors; 64 65 PetscFunctionBegin; 66 for (i = 0; i < n; i++) PetscCall((*ls->monitorftns[i])(ls, ls->monitorcontext[i])); 67 PetscFunctionReturn(PETSC_SUCCESS); 68 } 69 70 /*@C 71 SNESLineSearchMonitorSet - Sets an ADDITIONAL function that is to be used at every 72 iteration of the nonlinear solver to display the iteration's 73 progress. 74 75 Logically Collective 76 77 Input Parameters: 78 + ls - the `SNESLineSearch` context 79 . f - the monitor function 80 . mctx - [optional] user-defined context for private data for the 81 monitor routine (use `NULL` if no context is desired) 82 - monitordestroy - [optional] routine that frees monitor context 83 (may be `NULL`) 84 85 Level: intermediate 86 87 Note: 88 Several different monitoring routines may be set by calling 89 `SNESLineSearchMonitorSet()` multiple times; all will be called in the 90 order in which they were set. 91 92 Fortran Notes: 93 Only a single monitor function can be set for each `SNESLineSearch` object 94 95 .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchMonitorDefault()`, `SNESLineSearchMonitorCancel()` 96 @*/ 97 PetscErrorCode SNESLineSearchMonitorSet(SNESLineSearch ls, PetscErrorCode (*f)(SNESLineSearch, void *), void *mctx, PetscErrorCode (*monitordestroy)(void **)) 98 { 99 PetscInt i; 100 PetscBool identical; 101 102 PetscFunctionBegin; 103 PetscValidHeaderSpecific(ls, SNESLINESEARCH_CLASSID, 1); 104 for (i = 0; i < ls->numbermonitors; i++) { 105 PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))f, mctx, monitordestroy, (PetscErrorCode(*)(void))ls->monitorftns[i], ls->monitorcontext[i], ls->monitordestroy[i], &identical)); 106 if (identical) PetscFunctionReturn(PETSC_SUCCESS); 107 } 108 PetscCheck(ls->numbermonitors < MAXSNESLSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many monitors set"); 109 ls->monitorftns[ls->numbermonitors] = f; 110 ls->monitordestroy[ls->numbermonitors] = monitordestroy; 111 ls->monitorcontext[ls->numbermonitors++] = (void *)mctx; 112 PetscFunctionReturn(PETSC_SUCCESS); 113 } 114 115 /*@C 116 SNESLineSearchMonitorSolutionUpdate - Monitors each update of the function value the linesearch tries 117 118 Collective 119 120 Input Parameters: 121 + ls - the `SNES` linesearch object 122 - vf - the context for the monitor, in this case it is an `PetscViewerAndFormat` 123 124 Level: intermediate 125 126 .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESMonitorSet()`, `SNESMonitorSolution()` 127 @*/ 128 PetscErrorCode SNESLineSearchMonitorSolutionUpdate(SNESLineSearch ls, PetscViewerAndFormat *vf) 129 { 130 PetscViewer viewer = vf->viewer; 131 Vec Y, W, G; 132 133 PetscFunctionBegin; 134 PetscCall(SNESLineSearchGetVecs(ls, NULL, NULL, &Y, &W, &G)); 135 PetscCall(PetscViewerPushFormat(viewer, vf->format)); 136 PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted update to solution \n")); 137 PetscCall(VecView(Y, viewer)); 138 PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted new solution \n")); 139 PetscCall(VecView(W, viewer)); 140 PetscCall(PetscViewerASCIIPrintf(viewer, "LineSearch attempted updated function value\n")); 141 PetscCall(VecView(G, viewer)); 142 PetscCall(PetscViewerPopFormat(viewer)); 143 PetscFunctionReturn(PETSC_SUCCESS); 144 } 145 146 /*@ 147 SNESLineSearchCreate - Creates the line search context. 148 149 Logically Collective 150 151 Input Parameter: 152 . comm - MPI communicator for the line search (typically from the associated `SNES` context). 153 154 Output Parameter: 155 . outlinesearch - the new linesearch context 156 157 Level: developer 158 159 Note: 160 The preferred calling sequence for users is to use `SNESGetLineSearch()` to acquire the `SNESLineSearch` instance 161 already associated with the `SNES`. 162 163 .seealso: `SNESLineSearch`, `LineSearchDestroy()`, `SNESGetLineSearch()` 164 @*/ 165 PetscErrorCode SNESLineSearchCreate(MPI_Comm comm, SNESLineSearch *outlinesearch) 166 { 167 SNESLineSearch linesearch; 168 169 PetscFunctionBegin; 170 PetscAssertPointer(outlinesearch, 2); 171 PetscCall(SNESInitializePackage()); 172 *outlinesearch = NULL; 173 174 PetscCall(PetscHeaderCreate(linesearch, SNESLINESEARCH_CLASSID, "SNESLineSearch", "Linesearch", "SNESLineSearch", comm, SNESLineSearchDestroy, SNESLineSearchView)); 175 176 linesearch->vec_sol_new = NULL; 177 linesearch->vec_func_new = NULL; 178 linesearch->vec_sol = NULL; 179 linesearch->vec_func = NULL; 180 linesearch->vec_update = NULL; 181 182 linesearch->lambda = 1.0; 183 linesearch->fnorm = 1.0; 184 linesearch->ynorm = 1.0; 185 linesearch->xnorm = 1.0; 186 linesearch->result = SNES_LINESEARCH_SUCCEEDED; 187 linesearch->norms = PETSC_TRUE; 188 linesearch->keeplambda = PETSC_FALSE; 189 linesearch->damping = 1.0; 190 linesearch->maxstep = 1e8; 191 linesearch->steptol = 1e-12; 192 linesearch->rtol = 1e-8; 193 linesearch->atol = 1e-15; 194 linesearch->ltol = 1e-8; 195 linesearch->precheckctx = NULL; 196 linesearch->postcheckctx = NULL; 197 linesearch->max_its = 1; 198 linesearch->setupcalled = PETSC_FALSE; 199 linesearch->monitor = NULL; 200 *outlinesearch = linesearch; 201 PetscFunctionReturn(PETSC_SUCCESS); 202 } 203 204 /*@ 205 SNESLineSearchSetUp - Prepares the line search for being applied by allocating 206 any required vectors. 207 208 Collective 209 210 Input Parameter: 211 . linesearch - The `SNESLineSearch` instance. 212 213 Level: advanced 214 215 Note: 216 For most cases, this needn't be called by users or outside of `SNESLineSearchApply()`. 217 The only current case where this is called outside of this is for the VI 218 solvers, which modify the solution and work vectors before the first call 219 of `SNESLineSearchApply()`, requiring the `SNESLineSearch` work vectors to be 220 allocated upfront. 221 222 .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchReset()` 223 @*/ 224 PetscErrorCode SNESLineSearchSetUp(SNESLineSearch linesearch) 225 { 226 PetscFunctionBegin; 227 if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC)); 228 if (!linesearch->setupcalled) { 229 if (!linesearch->vec_sol_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_sol_new)); 230 if (!linesearch->vec_func_new) PetscCall(VecDuplicate(linesearch->vec_sol, &linesearch->vec_func_new)); 231 if (linesearch->ops->setup) PetscUseTypeMethod(linesearch, setup); 232 if (!linesearch->ops->snesfunc) PetscCall(SNESLineSearchSetFunction(linesearch, SNESComputeFunction)); 233 linesearch->lambda = linesearch->damping; 234 linesearch->setupcalled = PETSC_TRUE; 235 } 236 PetscFunctionReturn(PETSC_SUCCESS); 237 } 238 239 /*@ 240 SNESLineSearchReset - Undoes the `SNESLineSearchSetUp()` and deletes any `Vec`s or `Mat`s allocated by the line search. 241 242 Collective 243 244 Input Parameter: 245 . linesearch - The `SNESLineSearch` instance. 246 247 Level: developer 248 249 Note: 250 Usually only called by `SNESReset()` 251 252 .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetUp()` 253 @*/ 254 PetscErrorCode SNESLineSearchReset(SNESLineSearch linesearch) 255 { 256 PetscFunctionBegin; 257 if (linesearch->ops->reset) PetscUseTypeMethod(linesearch, reset); 258 259 PetscCall(VecDestroy(&linesearch->vec_sol_new)); 260 PetscCall(VecDestroy(&linesearch->vec_func_new)); 261 262 PetscCall(VecDestroyVecs(linesearch->nwork, &linesearch->work)); 263 264 linesearch->nwork = 0; 265 linesearch->setupcalled = PETSC_FALSE; 266 PetscFunctionReturn(PETSC_SUCCESS); 267 } 268 269 /*@C 270 SNESLineSearchSetFunction - Sets the function evaluation used by the `SNES` line search 271 ` 272 273 Input Parameters: 274 + linesearch - the `SNESLineSearch` context 275 - func - function evaluation routine, this is usually the function provided with `SNESSetFunction()` 276 277 Level: developer 278 279 .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESSetFunction()` 280 @*/ 281 PetscErrorCode SNESLineSearchSetFunction(SNESLineSearch linesearch, PetscErrorCode (*func)(SNES, Vec, Vec)) 282 { 283 PetscFunctionBegin; 284 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 285 linesearch->ops->snesfunc = func; 286 PetscFunctionReturn(PETSC_SUCCESS); 287 } 288 289 /*@C 290 SNESLineSearchSetPreCheck - Sets a user function that is called after the initial search direction has been computed but 291 before the line search routine has been applied. Allows the user to adjust the result of (usually a linear solve) that 292 determined the search direction. 293 294 Logically Collective 295 296 Input Parameters: 297 + linesearch - the `SNESLineSearch` context 298 . func - [optional] function evaluation routine, for the calling sequence see `SNESLineSearchPreCheck()` 299 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 300 301 Level: intermediate 302 303 Note: 304 Use `SNESLineSearchSetPostCheck()` to change the step after the line search. 305 search is complete. 306 307 Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed. 308 309 .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`, 310 `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()` 311 312 @*/ 313 PetscErrorCode SNESLineSearchSetPreCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void *ctx) 314 { 315 PetscFunctionBegin; 316 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 317 if (func) linesearch->ops->precheck = func; 318 if (ctx) linesearch->precheckctx = ctx; 319 PetscFunctionReturn(PETSC_SUCCESS); 320 } 321 322 /*@C 323 SNESLineSearchGetPreCheck - Gets the pre-check function for the line search routine. 324 325 Input Parameter: 326 . linesearch - the `SNESLineSearch` context 327 328 Output Parameters: 329 + func - [optional] function evaluation routine, for calling sequence see `SNESLineSearchPreCheck` 330 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 331 332 Level: intermediate 333 334 .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchSetPostCheck()` 335 @*/ 336 PetscErrorCode SNESLineSearchGetPreCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, PetscBool *, void *), void **ctx) 337 { 338 PetscFunctionBegin; 339 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 340 if (func) *func = linesearch->ops->precheck; 341 if (ctx) *ctx = linesearch->precheckctx; 342 PetscFunctionReturn(PETSC_SUCCESS); 343 } 344 345 /*@C 346 SNESLineSearchSetPostCheck - Sets a user function that is called after the line search has been applied to determine the step 347 direction and length. Allows the user a chance to change or override the decision of the line search routine 348 349 Logically Collective 350 351 Input Parameters: 352 + linesearch - the `SNESLineSearch` context 353 . func - [optional] function evaluation routine, for the calling sequence see `SNESLineSearchPostCheck` 354 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 355 356 Level: intermediate 357 358 Notes: 359 Use `SNESLineSearchSetPreCheck()` to change the step before the line search. 360 search is complete. 361 362 Use `SNESVISetVariableBounds()` and `SNESVISetComputeVariableBounds()` to cause `SNES` to automatically control the ranges of variables allowed. 363 364 .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchGetPostCheck()`, 365 `SNESVISetVariableBounds()`, `SNESVISetComputeVariableBounds()`, `SNESSetFunctionDomainError()`, `SNESSetJacobianDomainError()` 366 @*/ 367 PetscErrorCode SNESLineSearchSetPostCheck(SNESLineSearch linesearch, PetscErrorCode (*func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) 368 { 369 PetscFunctionBegin; 370 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 371 if (func) linesearch->ops->postcheck = func; 372 if (ctx) linesearch->postcheckctx = ctx; 373 PetscFunctionReturn(PETSC_SUCCESS); 374 } 375 376 /*@C 377 SNESLineSearchGetPostCheck - Gets the post-check function for the line search routine. 378 379 Input Parameter: 380 . linesearch - the `SNESLineSearch` context 381 382 Output Parameters: 383 + func - [optional] function evaluation routine, see for the calling sequence `SNESLineSearchPostCheck` 384 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 385 386 Level: intermediate 387 388 .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()` 389 @*/ 390 PetscErrorCode SNESLineSearchGetPostCheck(SNESLineSearch linesearch, PetscErrorCode (**func)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 391 { 392 PetscFunctionBegin; 393 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 394 if (func) *func = linesearch->ops->postcheck; 395 if (ctx) *ctx = linesearch->postcheckctx; 396 PetscFunctionReturn(PETSC_SUCCESS); 397 } 398 399 /*@ 400 SNESLineSearchPreCheck - Prepares the line search for being applied. 401 402 Logically Collective 403 404 Input Parameters: 405 + linesearch - The linesearch instance. 406 . X - The current solution 407 - Y - The step direction 408 409 Output Parameter: 410 . changed - Indicator that the precheck routine has changed anything 411 412 Level: advanced 413 414 Note: 415 This calls any function provided with `SNESLineSearchSetPreCheck()` 416 417 .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchPostCheck()`, `SNESLineSearchSetPreCheck()`, `SNESLineSearchGetPreCheck()`, `SNESLineSearchSetPostCheck()`, 418 `SNESLineSearchGetPostCheck()`` 419 @*/ 420 PetscErrorCode SNESLineSearchPreCheck(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed) 421 { 422 PetscFunctionBegin; 423 *changed = PETSC_FALSE; 424 if (linesearch->ops->precheck) { 425 PetscUseTypeMethod(linesearch, precheck, X, Y, changed, linesearch->precheckctx); 426 PetscValidLogicalCollectiveBool(linesearch, *changed, 4); 427 } 428 PetscFunctionReturn(PETSC_SUCCESS); 429 } 430 431 /*@ 432 SNESLineSearchPostCheck - Hook to modify step direction or updated solution after a successful linesearch 433 434 Logically Collective 435 436 Input Parameters: 437 + linesearch - The linesearch context 438 . X - The last solution 439 . Y - The step direction 440 - W - The updated solution, W = X + lambda*Y for some lambda 441 442 Output Parameters: 443 + changed_Y - Indicator if the direction Y has been changed. 444 - changed_W - Indicator if the new candidate solution W has been changed. 445 446 Level: developer 447 448 Note: 449 This calls any function provided with `SNESLineSearchSetPreCheck()` 450 451 .seealso: `SNESGetLineSearch()`, `SNESLineSearchPreCheck()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchSetPrecheck()`, `SNESLineSearchGetPrecheck()` 452 @*/ 453 PetscErrorCode SNESLineSearchPostCheck(SNESLineSearch linesearch, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 454 { 455 PetscFunctionBegin; 456 *changed_Y = PETSC_FALSE; 457 *changed_W = PETSC_FALSE; 458 if (linesearch->ops->postcheck) { 459 PetscUseTypeMethod(linesearch, postcheck, X, Y, W, changed_Y, changed_W, linesearch->postcheckctx); 460 PetscValidLogicalCollectiveBool(linesearch, *changed_Y, 5); 461 PetscValidLogicalCollectiveBool(linesearch, *changed_W, 6); 462 } 463 PetscFunctionReturn(PETSC_SUCCESS); 464 } 465 466 /*@C 467 SNESLineSearchPreCheckPicard - Implements a correction that is sometimes useful to improve the convergence rate of Picard iteration 468 469 Logically Collective 470 471 Input Parameters: 472 + linesearch - linesearch context 473 . X - base state for this step 474 - ctx - context for this function 475 476 Input/Output Parameter: 477 . Y - correction, possibly modified 478 479 Output Parameter: 480 . changed - flag indicating that Y was modified 481 482 Options Database Key: 483 + -snes_linesearch_precheck_picard - activate this routine 484 - -snes_linesearch_precheck_picard_angle - angle 485 486 Level: advanced 487 488 Notes: 489 This function should be passed to `SNESLineSearchSetPreCheck()` 490 491 The justification for this method involves the linear convergence of a Picard iteration 492 so the Picard linearization should be provided in place of the "Jacobian". This correction 493 is generally not useful when using a Newton linearization. 494 495 References: 496 . - * - Hindmarsh and Payne (1996) Time step limits for stable solutions of the ice sheet equation, Annals of Glaciology. 497 498 .seealso: `SNESLineSearch`, `SNESSetPicard()`, `SNESGetLineSearch()`, `SNESLineSearchSetPreCheck()` 499 @*/ 500 PetscErrorCode SNESLineSearchPreCheckPicard(SNESLineSearch linesearch, Vec X, Vec Y, PetscBool *changed, void *ctx) 501 { 502 PetscReal angle = *(PetscReal *)linesearch->precheckctx; 503 Vec Ylast; 504 PetscScalar dot; 505 PetscInt iter; 506 PetscReal ynorm, ylastnorm, theta, angle_radians; 507 SNES snes; 508 509 PetscFunctionBegin; 510 PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 511 PetscCall(PetscObjectQuery((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject *)&Ylast)); 512 if (!Ylast) { 513 PetscCall(VecDuplicate(Y, &Ylast)); 514 PetscCall(PetscObjectCompose((PetscObject)snes, "SNESLineSearchPreCheckPicard_Ylast", (PetscObject)Ylast)); 515 PetscCall(PetscObjectDereference((PetscObject)Ylast)); 516 } 517 PetscCall(SNESGetIterationNumber(snes, &iter)); 518 if (iter < 2) { 519 PetscCall(VecCopy(Y, Ylast)); 520 *changed = PETSC_FALSE; 521 PetscFunctionReturn(PETSC_SUCCESS); 522 } 523 524 PetscCall(VecDot(Y, Ylast, &dot)); 525 PetscCall(VecNorm(Y, NORM_2, &ynorm)); 526 PetscCall(VecNorm(Ylast, NORM_2, &ylastnorm)); 527 /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */ 528 theta = PetscAcosReal((PetscReal)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm), -1.0, 1.0)); 529 angle_radians = angle * PETSC_PI / 180.; 530 if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) { 531 /* Modify the step Y */ 532 PetscReal alpha, ydiffnorm; 533 PetscCall(VecAXPY(Ylast, -1.0, Y)); 534 PetscCall(VecNorm(Ylast, NORM_2, &ydiffnorm)); 535 alpha = (ydiffnorm > .001 * ylastnorm) ? ylastnorm / ydiffnorm : 1000.0; 536 PetscCall(VecCopy(Y, Ylast)); 537 PetscCall(VecScale(Y, alpha)); 538 PetscCall(PetscInfo(snes, "Angle %14.12e degrees less than threshold %14.12e, corrected step by alpha=%14.12e\n", (double)(theta * 180. / PETSC_PI), (double)angle, (double)alpha)); 539 *changed = PETSC_TRUE; 540 } else { 541 PetscCall(PetscInfo(snes, "Angle %14.12e degrees exceeds threshold %14.12e, no correction applied\n", (double)(theta * 180. / PETSC_PI), (double)angle)); 542 PetscCall(VecCopy(Y, Ylast)); 543 *changed = PETSC_FALSE; 544 } 545 PetscFunctionReturn(PETSC_SUCCESS); 546 } 547 548 /*@ 549 SNESLineSearchApply - Computes the line-search update. 550 551 Collective 552 553 Input Parameters: 554 + linesearch - The linesearch context 555 - Y - The search direction 556 557 Input/Output Parameters: 558 + X - The current solution, on output the new solution 559 . F - The current function, on output the new function 560 - fnorm - The current norm, on output the new function norm 561 562 Options Database Keys: 563 + -snes_linesearch_type - basic (or equivalently none), bt, l2, cp, nleqerr, shell 564 . -snes_linesearch_monitor [:filename] - Print progress of line searches 565 . -snes_linesearch_damping - The linesearch damping parameter, default is 1.0 (no damping) 566 . -snes_linesearch_norms - Turn on/off the linesearch norms computation (SNESLineSearchSetComputeNorms()) 567 . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess 568 - -snes_linesearch_max_it - The number of iterations for iterative line searches 569 570 Level: intermediate 571 572 Notes: 573 This is typically called from within a `SNESSolve()` implementation in order to 574 help with convergence of the nonlinear method. Various `SNES` types use line searches 575 in different ways, but the overarching theme is that a line search is used to determine 576 an optimal damping parameter of a step at each iteration of the method. Each 577 application of the line search may invoke `SNESComputeFunction()` several times, and 578 therefore may be fairly expensive. 579 580 .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchPreCheck()`, `SNESLineSearchPostCheck()`, `SNESSolve()`, `SNESComputeFunction()`, `SNESLineSearchSetComputeNorms()`, 581 `SNESLineSearchType`, `SNESLineSearchSetType()` 582 @*/ 583 PetscErrorCode SNESLineSearchApply(SNESLineSearch linesearch, Vec X, Vec F, PetscReal *fnorm, Vec Y) 584 { 585 PetscFunctionBegin; 586 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 587 PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 588 PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 589 PetscValidHeaderSpecific(Y, VEC_CLASSID, 5); 590 591 linesearch->result = SNES_LINESEARCH_SUCCEEDED; 592 593 linesearch->vec_sol = X; 594 linesearch->vec_update = Y; 595 linesearch->vec_func = F; 596 597 PetscCall(SNESLineSearchSetUp(linesearch)); 598 599 if (!linesearch->keeplambda) linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */ 600 601 if (fnorm) linesearch->fnorm = *fnorm; 602 else PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm)); 603 604 PetscCall(PetscLogEventBegin(SNESLINESEARCH_Apply, linesearch, X, F, Y)); 605 606 PetscUseTypeMethod(linesearch, apply); 607 608 PetscCall(PetscLogEventEnd(SNESLINESEARCH_Apply, linesearch, X, F, Y)); 609 610 if (fnorm) *fnorm = linesearch->fnorm; 611 PetscFunctionReturn(PETSC_SUCCESS); 612 } 613 614 /*@ 615 SNESLineSearchDestroy - Destroys the line search instance. 616 617 Collective 618 619 Input Parameter: 620 . linesearch - The linesearch context 621 622 Level: developer 623 624 .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchReset()`, `SNESDestroy()` 625 @*/ 626 PetscErrorCode SNESLineSearchDestroy(SNESLineSearch *linesearch) 627 { 628 PetscFunctionBegin; 629 if (!*linesearch) PetscFunctionReturn(PETSC_SUCCESS); 630 PetscValidHeaderSpecific((*linesearch), SNESLINESEARCH_CLASSID, 1); 631 if (--((PetscObject)(*linesearch))->refct > 0) { 632 *linesearch = NULL; 633 PetscFunctionReturn(PETSC_SUCCESS); 634 } 635 PetscCall(PetscObjectSAWsViewOff((PetscObject)*linesearch)); 636 PetscCall(SNESLineSearchReset(*linesearch)); 637 PetscTryTypeMethod(*linesearch, destroy); 638 PetscCall(PetscViewerDestroy(&(*linesearch)->monitor)); 639 PetscCall(SNESLineSearchMonitorCancel((*linesearch))); 640 PetscCall(PetscHeaderDestroy(linesearch)); 641 PetscFunctionReturn(PETSC_SUCCESS); 642 } 643 644 /*@ 645 SNESLineSearchSetDefaultMonitor - Turns on/off printing useful information and debugging output about the line search. 646 647 Logically Collective 648 649 Input Parameters: 650 + linesearch - the linesearch object 651 - viewer - an `PETSCVIEWERASCII` `PetscViewer` or `NULL` to turn off monitor 652 653 Options Database Key: 654 . -snes_linesearch_monitor [:filename] - enables the monitor 655 656 Level: intermediate 657 658 Developer Notes: 659 This monitor is implemented differently than the other line search monitors that are set with 660 `SNESLineSearchMonitorSet()` since it is called in many locations of the line search routines to display aspects of the 661 line search that are not visible to the other monitors. 662 663 .seealso: `SNESLineSearch`, `PETSCVIEWERASCII`, `SNESGetLineSearch()`, `SNESLineSearchGetDefaultMonitor()`, `PetscViewer`, `SNESLineSearchSetMonitor()`, 664 `SNESLineSearchMonitorSetFromOptions()` 665 @*/ 666 PetscErrorCode SNESLineSearchSetDefaultMonitor(SNESLineSearch linesearch, PetscViewer viewer) 667 { 668 PetscFunctionBegin; 669 if (viewer) PetscCall(PetscObjectReference((PetscObject)viewer)); 670 PetscCall(PetscViewerDestroy(&linesearch->monitor)); 671 linesearch->monitor = viewer; 672 PetscFunctionReturn(PETSC_SUCCESS); 673 } 674 675 /*@ 676 SNESLineSearchGetDefaultMonitor - Gets the `PetscViewer` instance for the line search monitor. 677 678 Logically Collective 679 680 Input Parameter: 681 . linesearch - linesearch context 682 683 Output Parameter: 684 . monitor - monitor context 685 686 Level: intermediate 687 688 .seealso: `SNESLineSearch`, `SNESGetLineSearch()`, `SNESLineSearchSetDefaultMonitor()`, `PetscViewer` 689 @*/ 690 PetscErrorCode SNESLineSearchGetDefaultMonitor(SNESLineSearch linesearch, PetscViewer *monitor) 691 { 692 PetscFunctionBegin; 693 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 694 *monitor = linesearch->monitor; 695 PetscFunctionReturn(PETSC_SUCCESS); 696 } 697 698 /*@C 699 SNESLineSearchMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 700 701 Collective 702 703 Input Parameters: 704 + ls - `SNESLineSearch` object you wish to monitor 705 . name - the monitor type one is seeking 706 . help - message indicating what monitoring is done 707 . manual - manual page for the monitor 708 . monitor - the monitor function 709 - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the `SNESLineSearch` or `PetscViewer` 710 711 Level: developer 712 713 .seealso: `SNESLineSearch`, `SNESLineSearchSetMonitor()`, `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, 714 `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 715 `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, 716 `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 717 `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 718 `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 719 `PetscOptionsFList()`, `PetscOptionsEList()` 720 @*/ 721 PetscErrorCode SNESLineSearchMonitorSetFromOptions(SNESLineSearch ls, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(SNESLineSearch, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(SNESLineSearch, PetscViewerAndFormat *)) 722 { 723 PetscViewer viewer; 724 PetscViewerFormat format; 725 PetscBool flg; 726 727 PetscFunctionBegin; 728 PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ls), ((PetscObject)ls)->options, ((PetscObject)ls)->prefix, name, &viewer, &format, &flg)); 729 if (flg) { 730 PetscViewerAndFormat *vf; 731 PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf)); 732 PetscCall(PetscObjectDereference((PetscObject)viewer)); 733 if (monitorsetup) PetscCall((*monitorsetup)(ls, vf)); 734 PetscCall(SNESLineSearchMonitorSet(ls, (PetscErrorCode(*)(SNESLineSearch, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy)); 735 } 736 PetscFunctionReturn(PETSC_SUCCESS); 737 } 738 739 /*@ 740 SNESLineSearchSetFromOptions - Sets options for the line search 741 742 Logically Collective 743 744 Input Parameter: 745 . linesearch - linesearch context 746 747 Options Database Keys: 748 + -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell 749 . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (bt supports 2 or 3) 750 . -snes_linesearch_norms - Turn on/off the linesearch norms for the basic linesearch typem (`SNESLineSearchSetComputeNorms()`) 751 . -snes_linesearch_minlambda - The minimum step length 752 . -snes_linesearch_maxstep - The maximum step size 753 . -snes_linesearch_rtol - Relative tolerance for iterative line searches 754 . -snes_linesearch_atol - Absolute tolerance for iterative line searches 755 . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches 756 . -snes_linesearch_max_it - The number of iterations for iterative line searches 757 . -snes_linesearch_monitor [:filename] - Print progress of line searches 758 . -snes_linesearch_monitor_solution_update [viewer:filename:format] - view each update tried by line search routine 759 . -snes_linesearch_damping - The linesearch damping parameter 760 . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess. 761 . -snes_linesearch_precheck_picard - Use precheck that speeds up convergence of picard method 762 - -snes_linesearch_precheck_picard_angle - Angle used in Picard precheck method 763 764 Level: intermediate 765 766 .seealso: `SNESLineSearch`, `SNESLineSearchCreate()`, `SNESLineSearchSetOrder()`, `SNESLineSearchSetType()`, `SNESLineSearchSetTolerances()`, `SNESLineSearchSetDamping()`, `SNESLineSearchPreCheckPicard()`, 767 `SNESLineSearchType`, `SNESLineSearchSetComputeNorms()` 768 @*/ 769 PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) 770 { 771 const char *deft = SNESLINESEARCHBASIC; 772 char type[256]; 773 PetscBool flg, set; 774 PetscViewer viewer; 775 776 PetscFunctionBegin; 777 PetscCall(SNESLineSearchRegisterAll()); 778 779 PetscObjectOptionsBegin((PetscObject)linesearch); 780 if (((PetscObject)linesearch)->type_name) deft = ((PetscObject)linesearch)->type_name; 781 PetscCall(PetscOptionsFList("-snes_linesearch_type", "Linesearch type", "SNESLineSearchSetType", SNESLineSearchList, deft, type, 256, &flg)); 782 if (flg) { 783 PetscCall(SNESLineSearchSetType(linesearch, type)); 784 } else if (!((PetscObject)linesearch)->type_name) { 785 PetscCall(SNESLineSearchSetType(linesearch, deft)); 786 } 787 788 PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)linesearch), ((PetscObject)linesearch)->options, ((PetscObject)linesearch)->prefix, "-snes_linesearch_monitor", &viewer, NULL, &set)); 789 if (set) { 790 PetscCall(SNESLineSearchSetDefaultMonitor(linesearch, viewer)); 791 PetscCall(PetscViewerDestroy(&viewer)); 792 } 793 PetscCall(SNESLineSearchMonitorSetFromOptions(linesearch, "-snes_linesearch_monitor_solution_update", "View correction at each iteration", "SNESLineSearchMonitorSolutionUpdate", SNESLineSearchMonitorSolutionUpdate, NULL)); 794 795 /* tolerances */ 796 PetscCall(PetscOptionsReal("-snes_linesearch_minlambda", "Minimum step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL)); 797 PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, NULL)); 798 PetscCall(PetscOptionsReal("-snes_linesearch_rtol", "Relative tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->rtol, &linesearch->rtol, NULL)); 799 PetscCall(PetscOptionsReal("-snes_linesearch_atol", "Absolute tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->atol, &linesearch->atol, NULL)); 800 PetscCall(PetscOptionsReal("-snes_linesearch_ltol", "Change in lambda tolerance for iterative line search", "SNESLineSearchSetTolerances", linesearch->ltol, &linesearch->ltol, NULL)); 801 PetscCall(PetscOptionsInt("-snes_linesearch_max_it", "Maximum iterations for iterative line searches", "SNESLineSearchSetTolerances", linesearch->max_its, &linesearch->max_its, NULL)); 802 803 /* damping parameters */ 804 PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL)); 805 806 PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL)); 807 808 /* precheck */ 809 PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set)); 810 if (set) { 811 if (flg) { 812 linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 813 814 PetscCall(PetscOptionsReal("-snes_linesearch_precheck_picard_angle", "Maximum angle at which to activate the correction", "none", linesearch->precheck_picard_angle, &linesearch->precheck_picard_angle, NULL)); 815 PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle)); 816 } else { 817 PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL)); 818 } 819 } 820 PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL)); 821 PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL)); 822 823 PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject); 824 825 PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject)); 826 PetscOptionsEnd(); 827 PetscFunctionReturn(PETSC_SUCCESS); 828 } 829 830 /*@ 831 SNESLineSearchView - Prints useful information about the line search 832 833 Logically Collective 834 835 Input Parameters: 836 + linesearch - linesearch context 837 - viewer - the viewer to display the line search information 838 839 Level: intermediate 840 841 .seealso: `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()` 842 @*/ 843 PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) 844 { 845 PetscBool iascii; 846 847 PetscFunctionBegin; 848 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 849 if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer)); 850 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 851 PetscCheckSameComm(linesearch, 1, viewer, 2); 852 853 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 854 if (iascii) { 855 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer)); 856 PetscCall(PetscViewerASCIIPushTab(viewer)); 857 PetscTryTypeMethod(linesearch, view, viewer); 858 PetscCall(PetscViewerASCIIPopTab(viewer)); 859 PetscCall(PetscViewerASCIIPrintf(viewer, " maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol)); 860 PetscCall(PetscViewerASCIIPrintf(viewer, " tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol)); 861 PetscCall(PetscViewerASCIIPrintf(viewer, " maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its)); 862 if (linesearch->ops->precheck) { 863 if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) { 864 PetscCall(PetscViewerASCIIPrintf(viewer, " using precheck step to speed up Picard convergence\n")); 865 } else { 866 PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined precheck step\n")); 867 } 868 } 869 if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined postcheck step\n")); 870 } 871 PetscFunctionReturn(PETSC_SUCCESS); 872 } 873 874 /*@C 875 SNESLineSearchGetType - Gets the linesearch type 876 877 Logically Collective 878 879 Input Parameter: 880 . linesearch - linesearch context 881 882 Output Parameter: 883 . type - The type of line search, or `NULL` if not set 884 885 Level: intermediate 886 887 .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()` 888 @*/ 889 PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type) 890 { 891 PetscFunctionBegin; 892 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 893 PetscAssertPointer(type, 2); 894 *type = ((PetscObject)linesearch)->type_name; 895 PetscFunctionReturn(PETSC_SUCCESS); 896 } 897 898 /*@C 899 SNESLineSearchSetType - Sets the linesearch type 900 901 Logically Collective 902 903 Input Parameters: 904 + linesearch - linesearch context 905 - type - The type of line search to be used, see `SNESLineSearchType` 906 907 Options Database Key: 908 . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell 909 910 Level: intermediate 911 912 .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()` 913 @*/ 914 PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) 915 { 916 PetscBool match; 917 PetscErrorCode (*r)(SNESLineSearch); 918 919 PetscFunctionBegin; 920 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 921 PetscAssertPointer(type, 2); 922 923 PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match)); 924 if (match) PetscFunctionReturn(PETSC_SUCCESS); 925 926 PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r)); 927 PetscCheck(r, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type); 928 /* Destroy the previous private linesearch context */ 929 if (linesearch->ops->destroy) { 930 PetscCall((*(linesearch)->ops->destroy)(linesearch)); 931 linesearch->ops->destroy = NULL; 932 } 933 /* Reinitialize function pointers in SNESLineSearchOps structure */ 934 linesearch->ops->apply = NULL; 935 linesearch->ops->view = NULL; 936 linesearch->ops->setfromoptions = NULL; 937 linesearch->ops->destroy = NULL; 938 939 PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type)); 940 PetscCall((*r)(linesearch)); 941 PetscFunctionReturn(PETSC_SUCCESS); 942 } 943 944 /*@ 945 SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation. 946 947 Input Parameters: 948 + linesearch - linesearch context 949 - snes - The snes instance 950 951 Level: developer 952 953 Note: 954 This happens automatically when the line search is obtained/created with 955 `SNESGetLineSearch()`. This routine is therefore mainly called within `SNES` 956 implementations. 957 958 .seealso: `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()`, `SNES` 959 @*/ 960 PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) 961 { 962 PetscFunctionBegin; 963 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 964 PetscValidHeaderSpecific(snes, SNES_CLASSID, 2); 965 linesearch->snes = snes; 966 PetscFunctionReturn(PETSC_SUCCESS); 967 } 968 969 /*@ 970 SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search. 971 972 Not Collective 973 974 Input Parameter: 975 . linesearch - linesearch context 976 977 Output Parameter: 978 . snes - The `SNES` instance 979 980 Level: developer 981 982 .seealso: `SNESLineSearch`, `SNESType`, `SNESLineSearchSetVecs()`, `SNES` 983 @*/ 984 PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) 985 { 986 PetscFunctionBegin; 987 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 988 PetscAssertPointer(snes, 2); 989 *snes = linesearch->snes; 990 PetscFunctionReturn(PETSC_SUCCESS); 991 } 992 993 /*@ 994 SNESLineSearchGetLambda - Gets the last linesearch steplength discovered. 995 996 Not Collective 997 998 Input Parameter: 999 . linesearch - linesearch context 1000 1001 Output Parameter: 1002 . lambda - The last steplength computed during `SNESLineSearchApply()` 1003 1004 Level: advanced 1005 1006 Note: 1007 This is useful in methods where the solver is ill-scaled and 1008 requires some adaptive notion of the difference in scale between the 1009 solution and the function. For instance, `SNESQN` may be scaled by the 1010 line search lambda using the argument -snes_qn_scaling ls. 1011 1012 .seealso: `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()` 1013 @*/ 1014 PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda) 1015 { 1016 PetscFunctionBegin; 1017 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1018 PetscAssertPointer(lambda, 2); 1019 *lambda = linesearch->lambda; 1020 PetscFunctionReturn(PETSC_SUCCESS); 1021 } 1022 1023 /*@ 1024 SNESLineSearchSetLambda - Sets the linesearch steplength 1025 1026 Input Parameters: 1027 + linesearch - linesearch context 1028 - lambda - The last steplength. 1029 1030 Level: advanced 1031 1032 Note: 1033 This routine is typically used within implementations of `SNESLineSearchApply()` 1034 to set the final steplength. This routine (and `SNESLineSearchGetLambda()`) were 1035 added in order to facilitate Quasi-Newton methods that use the previous steplength 1036 as an inner scaling parameter. 1037 1038 .seealso: `SNESLineSearch`, `SNESLineSearchGetLambda()` 1039 @*/ 1040 PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 1041 { 1042 PetscFunctionBegin; 1043 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1044 linesearch->lambda = lambda; 1045 PetscFunctionReturn(PETSC_SUCCESS); 1046 } 1047 1048 /*@ 1049 SNESLineSearchGetTolerances - Gets the tolerances for the linesearch. 1050 1051 Not Collective 1052 1053 Input Parameter: 1054 . linesearch - linesearch context 1055 1056 Output Parameters: 1057 + steptol - The minimum steplength 1058 . maxstep - The maximum steplength 1059 . rtol - The relative tolerance for iterative line searches 1060 . atol - The absolute tolerance for iterative line searches 1061 . ltol - The change in lambda tolerance for iterative line searches 1062 - max_its - The maximum number of iterations of the line search 1063 1064 Level: intermediate 1065 1066 Note: 1067 Different line searches may implement these parameters slightly differently as 1068 the type requires. 1069 1070 .seealso: `SNESLineSearch`, `SNESLineSearchSetTolerances()` 1071 @*/ 1072 PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 1073 { 1074 PetscFunctionBegin; 1075 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1076 if (steptol) { 1077 PetscAssertPointer(steptol, 2); 1078 *steptol = linesearch->steptol; 1079 } 1080 if (maxstep) { 1081 PetscAssertPointer(maxstep, 3); 1082 *maxstep = linesearch->maxstep; 1083 } 1084 if (rtol) { 1085 PetscAssertPointer(rtol, 4); 1086 *rtol = linesearch->rtol; 1087 } 1088 if (atol) { 1089 PetscAssertPointer(atol, 5); 1090 *atol = linesearch->atol; 1091 } 1092 if (ltol) { 1093 PetscAssertPointer(ltol, 6); 1094 *ltol = linesearch->ltol; 1095 } 1096 if (max_its) { 1097 PetscAssertPointer(max_its, 7); 1098 *max_its = linesearch->max_its; 1099 } 1100 PetscFunctionReturn(PETSC_SUCCESS); 1101 } 1102 1103 /*@ 1104 SNESLineSearchSetTolerances - Sets the tolerances for the linesearch. 1105 1106 Collective 1107 1108 Input Parameters: 1109 + linesearch - linesearch context 1110 . steptol - The minimum steplength 1111 . maxstep - The maximum steplength 1112 . rtol - The relative tolerance for iterative line searches 1113 . atol - The absolute tolerance for iterative line searches 1114 . ltol - The change in lambda tolerance for iterative line searches 1115 - max_its - The maximum number of iterations of the line search 1116 1117 Level: intermediate 1118 1119 Note: 1120 The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in 1121 place of an argument. 1122 1123 .seealso: `SNESLineSearch`, `SNESLineSearchGetTolerances()` 1124 @*/ 1125 PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 1126 { 1127 PetscFunctionBegin; 1128 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1129 PetscValidLogicalCollectiveReal(linesearch, steptol, 2); 1130 PetscValidLogicalCollectiveReal(linesearch, maxstep, 3); 1131 PetscValidLogicalCollectiveReal(linesearch, rtol, 4); 1132 PetscValidLogicalCollectiveReal(linesearch, atol, 5); 1133 PetscValidLogicalCollectiveReal(linesearch, ltol, 6); 1134 PetscValidLogicalCollectiveInt(linesearch, max_its, 7); 1135 1136 if (steptol != (PetscReal)PETSC_DEFAULT) { 1137 PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol); 1138 linesearch->steptol = steptol; 1139 } 1140 1141 if (maxstep != (PetscReal)PETSC_DEFAULT) { 1142 PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep); 1143 linesearch->maxstep = maxstep; 1144 } 1145 1146 if (rtol != (PetscReal)PETSC_DEFAULT) { 1147 PetscCheck(rtol >= 0.0 && rtol < 1.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Relative tolerance %14.12e must be non-negative and less than 1.0", (double)rtol); 1148 linesearch->rtol = rtol; 1149 } 1150 1151 if (atol != (PetscReal)PETSC_DEFAULT) { 1152 PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol); 1153 linesearch->atol = atol; 1154 } 1155 1156 if (ltol != (PetscReal)PETSC_DEFAULT) { 1157 PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol); 1158 linesearch->ltol = ltol; 1159 } 1160 1161 if (max_its != PETSC_DEFAULT) { 1162 PetscCheck(max_its >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_its); 1163 linesearch->max_its = max_its; 1164 } 1165 PetscFunctionReturn(PETSC_SUCCESS); 1166 } 1167 1168 /*@ 1169 SNESLineSearchGetDamping - Gets the line search damping parameter. 1170 1171 Input Parameter: 1172 . linesearch - linesearch context 1173 1174 Output Parameter: 1175 . damping - The damping parameter 1176 1177 Level: advanced 1178 1179 .seealso: `SNESLineSearchGetStepTolerance()`, `SNESQN` 1180 @*/ 1181 PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping) 1182 { 1183 PetscFunctionBegin; 1184 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1185 PetscAssertPointer(damping, 2); 1186 *damping = linesearch->damping; 1187 PetscFunctionReturn(PETSC_SUCCESS); 1188 } 1189 1190 /*@ 1191 SNESLineSearchSetDamping - Sets the line search damping parameter. 1192 1193 Input Parameters: 1194 + linesearch - linesearch context 1195 - damping - The damping parameter 1196 1197 Options Database Key: 1198 . -snes_linesearch_damping <damping> - the damping value 1199 1200 Level: intermediate 1201 1202 Note: 1203 The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter. 1204 The use of the damping parameter in the l2 and cp line searches is much more subtle; 1205 it is used as a starting point in calculating the secant step. However, the eventual 1206 step may be of greater length than the damping parameter. In the bt line search it is 1207 used as the maximum possible step length, as the bt line search only backtracks. 1208 1209 .seealso: `SNESLineSearch`, `SNESLineSearchGetDamping()` 1210 @*/ 1211 PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping) 1212 { 1213 PetscFunctionBegin; 1214 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1215 linesearch->damping = damping; 1216 PetscFunctionReturn(PETSC_SUCCESS); 1217 } 1218 1219 /*@ 1220 SNESLineSearchGetOrder - Gets the line search approximation order. 1221 1222 Input Parameter: 1223 . linesearch - linesearch context 1224 1225 Output Parameter: 1226 . order - The order 1227 1228 Level: intermediate 1229 1230 .seealso: `SNESLineSearch`, `SNESLineSearchSetOrder()` 1231 @*/ 1232 PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order) 1233 { 1234 PetscFunctionBegin; 1235 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1236 PetscAssertPointer(order, 2); 1237 *order = linesearch->order; 1238 PetscFunctionReturn(PETSC_SUCCESS); 1239 } 1240 1241 /*@ 1242 SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search 1243 1244 Input Parameters: 1245 + linesearch - linesearch context 1246 - order - The order 1247 1248 Level: intermediate 1249 1250 Values for `order`\: 1251 + 1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order 1252 . 2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order 1253 - 3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order 1254 1255 Note: 1256 These orders are supported by `SNESLINESEARCHBT` and `SNESLINESEARCHCP` 1257 1258 .seealso: `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()` 1259 @*/ 1260 PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order) 1261 { 1262 PetscFunctionBegin; 1263 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1264 linesearch->order = order; 1265 PetscFunctionReturn(PETSC_SUCCESS); 1266 } 1267 1268 /*@ 1269 SNESLineSearchGetNorms - Gets the norms for for X, Y, and F. 1270 1271 Not Collective 1272 1273 Input Parameter: 1274 . linesearch - linesearch context 1275 1276 Output Parameters: 1277 + xnorm - The norm of the current solution 1278 . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution. 1279 - ynorm - The norm of the current update 1280 1281 Level: developer 1282 1283 Note: 1284 Some values may not be up to date at particular points in the code. 1285 1286 This, in combination with `SNESLineSearchSetNorms()`, allow the line search and the `SNESSolve_XXX()` to share 1287 computed values. 1288 1289 .seealso: `SNESLineSearch`, `SNESLineSearchSetNorms()` `SNESLineSearchGetVecs()` 1290 @*/ 1291 PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm) 1292 { 1293 PetscFunctionBegin; 1294 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1295 if (xnorm) *xnorm = linesearch->xnorm; 1296 if (fnorm) *fnorm = linesearch->fnorm; 1297 if (ynorm) *ynorm = linesearch->ynorm; 1298 PetscFunctionReturn(PETSC_SUCCESS); 1299 } 1300 1301 /*@ 1302 SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 1303 1304 Collective 1305 1306 Input Parameters: 1307 + linesearch - linesearch context 1308 . xnorm - The norm of the current solution 1309 . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution 1310 - ynorm - The norm of the current update 1311 1312 Level: developer 1313 1314 .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1315 @*/ 1316 PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 1317 { 1318 PetscFunctionBegin; 1319 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1320 linesearch->xnorm = xnorm; 1321 linesearch->fnorm = fnorm; 1322 linesearch->ynorm = ynorm; 1323 PetscFunctionReturn(PETSC_SUCCESS); 1324 } 1325 1326 /*@ 1327 SNESLineSearchComputeNorms - Computes the norms of X, F, and Y. 1328 1329 Input Parameter: 1330 . linesearch - linesearch context 1331 1332 Options Database Key: 1333 . -snes_linesearch_norms - turn norm computation on or off 1334 1335 Level: intermediate 1336 1337 .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()` 1338 @*/ 1339 PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 1340 { 1341 SNES snes; 1342 1343 PetscFunctionBegin; 1344 if (linesearch->norms) { 1345 if (linesearch->ops->vinorm) { 1346 PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 1347 PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 1348 PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 1349 PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm)); 1350 } else { 1351 PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 1352 PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 1353 PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 1354 PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 1355 PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 1356 PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 1357 } 1358 } 1359 PetscFunctionReturn(PETSC_SUCCESS); 1360 } 1361 1362 /*@ 1363 SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 1364 1365 Input Parameters: 1366 + linesearch - linesearch context 1367 - flg - indicates whether or not to compute norms 1368 1369 Options Database Key: 1370 . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) linesearch 1371 1372 Level: intermediate 1373 1374 Note: 1375 This is most relevant to the `SNESLINESEARCHBASIC` (or equivalently `SNESLINESEARCHNONE`) line search type since most line searches have a stopping criteria involving the norm. 1376 1377 .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC` 1378 @*/ 1379 PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 1380 { 1381 PetscFunctionBegin; 1382 linesearch->norms = flg; 1383 PetscFunctionReturn(PETSC_SUCCESS); 1384 } 1385 1386 /*@ 1387 SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context 1388 1389 Not Collective on but the vectors are parallel 1390 1391 Input Parameter: 1392 . linesearch - linesearch context 1393 1394 Output Parameters: 1395 + X - Solution vector 1396 . F - Function vector 1397 . Y - Search direction vector 1398 . W - Solution work vector 1399 - G - Function work vector 1400 1401 Level: advanced 1402 1403 Notes: 1404 At the beginning of a line search application, `X` should contain a 1405 solution and the vector `F` the function computed at `X`. At the end of the 1406 line search application, `X` should contain the new solution, and `F` the 1407 function evaluated at the new solution. 1408 1409 These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller 1410 1411 .seealso: `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1412 @*/ 1413 PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G) 1414 { 1415 PetscFunctionBegin; 1416 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1417 if (X) { 1418 PetscAssertPointer(X, 2); 1419 *X = linesearch->vec_sol; 1420 } 1421 if (F) { 1422 PetscAssertPointer(F, 3); 1423 *F = linesearch->vec_func; 1424 } 1425 if (Y) { 1426 PetscAssertPointer(Y, 4); 1427 *Y = linesearch->vec_update; 1428 } 1429 if (W) { 1430 PetscAssertPointer(W, 5); 1431 *W = linesearch->vec_sol_new; 1432 } 1433 if (G) { 1434 PetscAssertPointer(G, 6); 1435 *G = linesearch->vec_func_new; 1436 } 1437 PetscFunctionReturn(PETSC_SUCCESS); 1438 } 1439 1440 /*@ 1441 SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context 1442 1443 Logically Collective 1444 1445 Input Parameters: 1446 + linesearch - linesearch context 1447 . X - Solution vector 1448 . F - Function vector 1449 . Y - Search direction vector 1450 . W - Solution work vector 1451 - G - Function work vector 1452 1453 Level: advanced 1454 1455 .seealso: `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()` 1456 @*/ 1457 PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G) 1458 { 1459 PetscFunctionBegin; 1460 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1461 if (X) { 1462 PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 1463 linesearch->vec_sol = X; 1464 } 1465 if (F) { 1466 PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 1467 linesearch->vec_func = F; 1468 } 1469 if (Y) { 1470 PetscValidHeaderSpecific(Y, VEC_CLASSID, 4); 1471 linesearch->vec_update = Y; 1472 } 1473 if (W) { 1474 PetscValidHeaderSpecific(W, VEC_CLASSID, 5); 1475 linesearch->vec_sol_new = W; 1476 } 1477 if (G) { 1478 PetscValidHeaderSpecific(G, VEC_CLASSID, 6); 1479 linesearch->vec_func_new = G; 1480 } 1481 PetscFunctionReturn(PETSC_SUCCESS); 1482 } 1483 1484 /*@C 1485 SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1486 `SNESLineSearch` options in the database. 1487 1488 Logically Collective 1489 1490 Input Parameters: 1491 + linesearch - the `SNESLineSearch` context 1492 - prefix - the prefix to prepend to all option names 1493 1494 Level: advanced 1495 1496 Note: 1497 A hyphen (-) must NOT be given at the beginning of the prefix name. 1498 The first character of all runtime options is AUTOMATICALLY the hyphen. 1499 1500 .seealso: `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()` 1501 @*/ 1502 PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[]) 1503 { 1504 PetscFunctionBegin; 1505 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1506 PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix)); 1507 PetscFunctionReturn(PETSC_SUCCESS); 1508 } 1509 1510 /*@C 1511 SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 1512 SNESLineSearch options in the database. 1513 1514 Not Collective 1515 1516 Input Parameter: 1517 . linesearch - the `SNESLineSearch` context 1518 1519 Output Parameter: 1520 . prefix - pointer to the prefix string used 1521 1522 Level: advanced 1523 1524 Fortran Notes: 1525 The user should pass in a string 'prefix' of 1526 sufficient length to hold the prefix. 1527 1528 .seealso: `SNESLineSearch`, `SNESAppendOptionsPrefix()` 1529 @*/ 1530 PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[]) 1531 { 1532 PetscFunctionBegin; 1533 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1534 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix)); 1535 PetscFunctionReturn(PETSC_SUCCESS); 1536 } 1537 1538 /*@C 1539 SNESLineSearchSetWorkVecs - Sets work vectors for the line search. 1540 1541 Input Parameters: 1542 + linesearch - the `SNESLineSearch` context 1543 - nwork - the number of work vectors 1544 1545 Level: developer 1546 1547 .seealso: `SNESLineSearch`, `SNESSetWorkVecs()` 1548 @*/ 1549 PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) 1550 { 1551 PetscFunctionBegin; 1552 PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 1553 PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work)); 1554 PetscFunctionReturn(PETSC_SUCCESS); 1555 } 1556 1557 /*@ 1558 SNESLineSearchGetReason - Gets the success/failure status of the last line search application 1559 1560 Input Parameter: 1561 . linesearch - linesearch context 1562 1563 Output Parameter: 1564 . result - The success or failure status 1565 1566 Level: developer 1567 1568 Note: 1569 This is typically called after `SNESLineSearchApply()` in order to determine if the line-search failed 1570 (and set the SNES convergence accordingly). 1571 1572 .seealso: `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason` 1573 @*/ 1574 PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) 1575 { 1576 PetscFunctionBegin; 1577 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1578 PetscAssertPointer(result, 2); 1579 *result = linesearch->result; 1580 PetscFunctionReturn(PETSC_SUCCESS); 1581 } 1582 1583 /*@ 1584 SNESLineSearchSetReason - Sets the success/failure status of the last line search application 1585 1586 Input Parameters: 1587 + linesearch - linesearch context 1588 - result - The success or failure status 1589 1590 Level: developer 1591 1592 Note: 1593 This is typically called in a `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set 1594 the success or failure of the line search method. 1595 1596 .seealso: `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()` 1597 @*/ 1598 PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) 1599 { 1600 PetscFunctionBegin; 1601 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1602 linesearch->result = result; 1603 PetscFunctionReturn(PETSC_SUCCESS); 1604 } 1605 1606 // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation 1607 /*@C 1608 SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 1609 1610 Logically Collective 1611 1612 Input Parameters: 1613 + linesearch - the linesearch object 1614 . projectfunc - function for projecting the function to the bounds 1615 - normfunc - function for computing the norm of an active set 1616 1617 Calling sequence of `projectfunc`: 1618 $ PetscErrorCode projectfunc(SNES snes, Vec X) 1619 + snes - nonlinear context 1620 - X - current solution, store the projected solution here 1621 1622 Calling sequence of `normfunc`: 1623 $ PetscErrorCode normfunc(SNES snes, Vec X, Vec F, PetscScalar *fnorm) 1624 + snes - nonlinear context 1625 . X - current solution 1626 . F - current residual 1627 - fnorm - the computed norm 1628 1629 Level: advanced 1630 1631 Notes: 1632 The VI solvers require projection of the solution to the feasible set. `projectfunc` should implement this. 1633 1634 The VI solvers require special evaluation of the function norm such that the norm is only calculated 1635 on the inactive set. This should be implemented by `normfunc`. 1636 1637 .seealso: `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` 1638 @*/ 1639 PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 1640 { 1641 PetscFunctionBegin; 1642 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1643 if (projectfunc) linesearch->ops->viproject = projectfunc; 1644 if (normfunc) linesearch->ops->vinorm = normfunc; 1645 PetscFunctionReturn(PETSC_SUCCESS); 1646 } 1647 1648 /*@C 1649 SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 1650 1651 Not Collective 1652 1653 Input Parameter: 1654 . linesearch - the line search context, obtain with `SNESGetLineSearch()` 1655 1656 Output Parameters: 1657 + projectfunc - function for projecting the function to the bounds 1658 - normfunc - function for computing the norm of an active set 1659 1660 Level: advanced 1661 1662 .seealso: `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()` 1663 @*/ 1664 PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 1665 { 1666 PetscFunctionBegin; 1667 if (projectfunc) *projectfunc = linesearch->ops->viproject; 1668 if (normfunc) *normfunc = linesearch->ops->vinorm; 1669 PetscFunctionReturn(PETSC_SUCCESS); 1670 } 1671 1672 /*@C 1673 SNESLineSearchRegister - register a line search type 1674 1675 Input Parameters: 1676 + sname - name of the line search 1677 - function - the creation function for that type 1678 1679 Calling sequence of `function`: 1680 . ls - the linesearch context 1681 1682 Level: advanced 1683 1684 .seealso: `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()` 1685 @*/ 1686 PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch ls)) 1687 { 1688 PetscFunctionBegin; 1689 PetscCall(SNESInitializePackage()); 1690 PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function)); 1691 PetscFunctionReturn(PETSC_SUCCESS); 1692 } 1693