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