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