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->maxstep = 1e8; 196 linesearch->steptol = 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_its = 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, l2, 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 search length 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 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, l2, cp, nleqerr, bisection, shell 793 . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (bt supports 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 step length 796 . -snes_linesearch_maxstep - The maximum step size 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 search length 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`, `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 step length", "SNESLineSearchSetTolerances", linesearch->steptol, &linesearch->steptol, NULL)); 838 PetscCall(PetscOptionsReal("-snes_linesearch_maxstep", "Maximum step size", "SNESLineSearchSetTolerances", linesearch->maxstep, &linesearch->maxstep, 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_its, &linesearch->max_its, NULL)); 843 844 /* damping parameters */ 845 PetscCall(PetscOptionsReal("-snes_linesearch_damping", "Line search damping and initial step guess", "SNESLineSearchSetDamping", linesearch->damping, &linesearch->damping, NULL)); 846 847 PetscCall(PetscOptionsBool("-snes_linesearch_keeplambda", "Use previous lambda as damping", "SNESLineSearchSetKeepLambda", linesearch->keeplambda, &linesearch->keeplambda, NULL)); 848 849 /* precheck */ 850 PetscCall(PetscOptionsBool("-snes_linesearch_precheck_picard", "Use a correction that sometimes improves convergence of Picard iteration", "SNESLineSearchPreCheckPicard", flg, &flg, &set)); 851 if (set) { 852 if (flg) { 853 linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 854 855 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)); 856 PetscCall(SNESLineSearchSetPreCheck(linesearch, SNESLineSearchPreCheckPicard, &linesearch->precheck_picard_angle)); 857 } else { 858 PetscCall(SNESLineSearchSetPreCheck(linesearch, NULL, NULL)); 859 } 860 } 861 PetscCall(PetscOptionsInt("-snes_linesearch_order", "Order of approximation used in the line search", "SNESLineSearchSetOrder", linesearch->order, &linesearch->order, NULL)); 862 PetscCall(PetscOptionsBool("-snes_linesearch_norms", "Compute final norms in line search", "SNESLineSearchSetComputeNorms", linesearch->norms, &linesearch->norms, NULL)); 863 864 PetscTryTypeMethod(linesearch, setfromoptions, PetscOptionsObject); 865 866 PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)linesearch, PetscOptionsObject)); 867 PetscOptionsEnd(); 868 PetscFunctionReturn(PETSC_SUCCESS); 869 } 870 871 /*@ 872 SNESLineSearchView - Prints useful information about the line search 873 874 Logically Collective 875 876 Input Parameters: 877 + linesearch - line search context 878 - viewer - the `PetscViewer` to display the line search information to 879 880 Level: intermediate 881 882 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `PetscViewer`, `SNESLineSearchCreate()` 883 @*/ 884 PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) 885 { 886 PetscBool iascii; 887 888 PetscFunctionBegin; 889 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 890 if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)linesearch), &viewer)); 891 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 892 PetscCheckSameComm(linesearch, 1, viewer, 2); 893 894 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 895 if (iascii) { 896 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)linesearch, viewer)); 897 PetscCall(PetscViewerASCIIPushTab(viewer)); 898 PetscTryTypeMethod(linesearch, view, viewer); 899 PetscCall(PetscViewerASCIIPopTab(viewer)); 900 PetscCall(PetscViewerASCIIPrintf(viewer, " maxstep=%e, minlambda=%e\n", (double)linesearch->maxstep, (double)linesearch->steptol)); 901 PetscCall(PetscViewerASCIIPrintf(viewer, " tolerances: relative=%e, absolute=%e, lambda=%e\n", (double)linesearch->rtol, (double)linesearch->atol, (double)linesearch->ltol)); 902 PetscCall(PetscViewerASCIIPrintf(viewer, " maximum iterations=%" PetscInt_FMT "\n", linesearch->max_its)); 903 if (linesearch->ops->precheck) { 904 if (linesearch->ops->precheck == SNESLineSearchPreCheckPicard) { 905 PetscCall(PetscViewerASCIIPrintf(viewer, " using precheck step to speed up Picard convergence\n")); 906 } else { 907 PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined precheck step\n")); 908 } 909 } 910 if (linesearch->ops->postcheck) PetscCall(PetscViewerASCIIPrintf(viewer, " using user-defined postcheck step\n")); 911 } 912 PetscFunctionReturn(PETSC_SUCCESS); 913 } 914 915 /*@ 916 SNESLineSearchGetType - Gets the `SNESLinesearchType` of a `SNESLineSearch` 917 918 Logically Collective 919 920 Input Parameter: 921 . linesearch - the line search context 922 923 Output Parameter: 924 . type - The type of line search, or `NULL` if not set 925 926 Level: intermediate 927 928 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchSetType()` 929 @*/ 930 PetscErrorCode SNESLineSearchGetType(SNESLineSearch linesearch, SNESLineSearchType *type) 931 { 932 PetscFunctionBegin; 933 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 934 PetscAssertPointer(type, 2); 935 *type = ((PetscObject)linesearch)->type_name; 936 PetscFunctionReturn(PETSC_SUCCESS); 937 } 938 939 /*@ 940 SNESLineSearchSetType - Sets the `SNESLinesearchType` of a `SNESLineSearch` object to indicate the line search algorithm that should be used by a given `SNES` solver 941 942 Logically Collective 943 944 Input Parameters: 945 + linesearch - the line search context 946 - type - The type of line search to be used, see `SNESLineSearchType` 947 948 Options Database Key: 949 . -snes_linesearch_type <type> - basic (or equivalently none), bt, l2, cp, nleqerr, shell 950 951 Level: intermediate 952 953 Note: 954 The `SNESLineSearch` object is generally obtained with `SNESGetLineSearch()` 955 956 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetFromOptions()`, `SNESLineSearchGetType()`, 957 `SNESGetLineSearch()` 958 @*/ 959 PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, SNESLineSearchType type) 960 { 961 PetscBool match; 962 PetscErrorCode (*r)(SNESLineSearch); 963 964 PetscFunctionBegin; 965 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 966 PetscAssertPointer(type, 2); 967 968 PetscCall(PetscObjectTypeCompare((PetscObject)linesearch, type, &match)); 969 if (match) PetscFunctionReturn(PETSC_SUCCESS); 970 971 PetscCall(PetscFunctionListFind(SNESLineSearchList, type, &r)); 972 PetscCheck(r, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Line Search type %s", type); 973 /* Destroy the previous private line search context */ 974 PetscTryTypeMethod(linesearch, destroy); 975 linesearch->ops->destroy = NULL; 976 /* Reinitialize function pointers in SNESLineSearchOps structure */ 977 linesearch->ops->apply = NULL; 978 linesearch->ops->view = NULL; 979 linesearch->ops->setfromoptions = NULL; 980 linesearch->ops->destroy = NULL; 981 982 PetscCall(PetscObjectChangeTypeName((PetscObject)linesearch, type)); 983 PetscCall((*r)(linesearch)); 984 PetscFunctionReturn(PETSC_SUCCESS); 985 } 986 987 /*@ 988 SNESLineSearchSetSNES - Sets the `SNES` for the linesearch for function evaluation. 989 990 Input Parameters: 991 + linesearch - the line search context 992 - snes - The `SNES` instance 993 994 Level: developer 995 996 Note: 997 This happens automatically when the line search is obtained/created with 998 `SNESGetLineSearch()`. This routine is therefore mainly called within `SNES` 999 implementations. 1000 1001 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetSNES()`, `SNESLineSearchSetVecs()` 1002 @*/ 1003 PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) 1004 { 1005 PetscFunctionBegin; 1006 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1007 PetscValidHeaderSpecific(snes, SNES_CLASSID, 2); 1008 linesearch->snes = snes; 1009 PetscFunctionReturn(PETSC_SUCCESS); 1010 } 1011 1012 /*@ 1013 SNESLineSearchGetSNES - Gets the `SNES` instance associated with the line search. 1014 1015 Not Collective 1016 1017 Input Parameter: 1018 . linesearch - the line search context 1019 1020 Output Parameter: 1021 . snes - The `SNES` instance 1022 1023 Level: developer 1024 1025 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESType`, `SNESLineSearchSetVecs()` 1026 @*/ 1027 PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) 1028 { 1029 PetscFunctionBegin; 1030 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1031 PetscAssertPointer(snes, 2); 1032 *snes = linesearch->snes; 1033 PetscFunctionReturn(PETSC_SUCCESS); 1034 } 1035 1036 /*@ 1037 SNESLineSearchGetLambda - Gets the last line search steplength used 1038 1039 Not Collective 1040 1041 Input Parameter: 1042 . linesearch - the line search context 1043 1044 Output Parameter: 1045 . lambda - The last steplength computed during `SNESLineSearchApply()` 1046 1047 Level: advanced 1048 1049 Note: 1050 This is useful in methods where the solver is ill-scaled and 1051 requires some adaptive notion of the difference in scale between the 1052 solution and the function. For instance, `SNESQN` may be scaled by the 1053 line search lambda using the argument -snes_qn_scaling ls. 1054 1055 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetLambda()`, `SNESLineSearchGetDamping()`, `SNESLineSearchApply()` 1056 @*/ 1057 PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch, PetscReal *lambda) 1058 { 1059 PetscFunctionBegin; 1060 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1061 PetscAssertPointer(lambda, 2); 1062 *lambda = linesearch->lambda; 1063 PetscFunctionReturn(PETSC_SUCCESS); 1064 } 1065 1066 /*@ 1067 SNESLineSearchSetLambda - Sets the line search steplength 1068 1069 Input Parameters: 1070 + linesearch - line search context 1071 - lambda - The steplength to use 1072 1073 Level: advanced 1074 1075 Note: 1076 This routine is typically used within implementations of `SNESLineSearchApply()` 1077 to set the final steplength. This routine (and `SNESLineSearchGetLambda()`) were 1078 added in order to facilitate Quasi-Newton methods that use the previous steplength 1079 as an inner scaling parameter. 1080 1081 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetLambda()` 1082 @*/ 1083 PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 1084 { 1085 PetscFunctionBegin; 1086 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1087 linesearch->lambda = lambda; 1088 PetscFunctionReturn(PETSC_SUCCESS); 1089 } 1090 1091 /*@ 1092 SNESLineSearchGetTolerances - Gets the tolerances for the line search. 1093 1094 Not Collective 1095 1096 Input Parameter: 1097 . linesearch - the line search context 1098 1099 Output Parameters: 1100 + steptol - The minimum steplength 1101 . maxstep - The maximum steplength 1102 . rtol - The relative tolerance for iterative line searches 1103 . atol - The absolute tolerance for iterative line searches 1104 . ltol - The change in lambda tolerance for iterative line searches 1105 - max_its - The maximum number of iterations of the line search 1106 1107 Level: intermediate 1108 1109 Note: 1110 Different line searches may implement these parameters slightly differently as 1111 the type requires. 1112 1113 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetTolerances()` 1114 @*/ 1115 PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch, PetscReal *steptol, PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 1116 { 1117 PetscFunctionBegin; 1118 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1119 if (steptol) { 1120 PetscAssertPointer(steptol, 2); 1121 *steptol = linesearch->steptol; 1122 } 1123 if (maxstep) { 1124 PetscAssertPointer(maxstep, 3); 1125 *maxstep = linesearch->maxstep; 1126 } 1127 if (rtol) { 1128 PetscAssertPointer(rtol, 4); 1129 *rtol = linesearch->rtol; 1130 } 1131 if (atol) { 1132 PetscAssertPointer(atol, 5); 1133 *atol = linesearch->atol; 1134 } 1135 if (ltol) { 1136 PetscAssertPointer(ltol, 6); 1137 *ltol = linesearch->ltol; 1138 } 1139 if (max_its) { 1140 PetscAssertPointer(max_its, 7); 1141 *max_its = linesearch->max_its; 1142 } 1143 PetscFunctionReturn(PETSC_SUCCESS); 1144 } 1145 1146 /*@ 1147 SNESLineSearchSetTolerances - Sets the tolerances for the linesearch. 1148 1149 Collective 1150 1151 Input Parameters: 1152 + linesearch - the line search context 1153 . steptol - The minimum steplength 1154 . maxstep - The maximum steplength 1155 . rtol - The relative tolerance for iterative line searches 1156 . atol - The absolute tolerance for iterative line searches 1157 . ltol - The change in lambda tolerance for iterative line searches 1158 - max_it - The maximum number of iterations of the line search 1159 1160 Options Database Keys: 1161 + -snes_linesearch_minlambda - The minimum step length 1162 . -snes_linesearch_maxstep - The maximum step size 1163 . -snes_linesearch_rtol - Relative tolerance for iterative line searches 1164 . -snes_linesearch_atol - Absolute tolerance for iterative line searches 1165 . -snes_linesearch_ltol - Change in lambda tolerance for iterative line searches 1166 - -snes_linesearch_max_it - The number of iterations for iterative line searches 1167 1168 Level: intermediate 1169 1170 Note: 1171 The user may choose to not set any of the tolerances using `PETSC_DEFAULT` in place of an argument. 1172 1173 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetTolerances()` 1174 @*/ 1175 PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch, PetscReal steptol, PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_it) 1176 { 1177 PetscFunctionBegin; 1178 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1179 PetscValidLogicalCollectiveReal(linesearch, steptol, 2); 1180 PetscValidLogicalCollectiveReal(linesearch, maxstep, 3); 1181 PetscValidLogicalCollectiveReal(linesearch, rtol, 4); 1182 PetscValidLogicalCollectiveReal(linesearch, atol, 5); 1183 PetscValidLogicalCollectiveReal(linesearch, ltol, 6); 1184 PetscValidLogicalCollectiveInt(linesearch, max_it, 7); 1185 1186 if (steptol != (PetscReal)PETSC_DEFAULT) { 1187 PetscCheck(steptol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Minimum step length %14.12e must be non-negative", (double)steptol); 1188 linesearch->steptol = steptol; 1189 } 1190 1191 if (maxstep != (PetscReal)PETSC_DEFAULT) { 1192 PetscCheck(maxstep >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum step length %14.12e must be non-negative", (double)maxstep); 1193 linesearch->maxstep = maxstep; 1194 } 1195 1196 if (rtol != (PetscReal)PETSC_DEFAULT) { 1197 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); 1198 linesearch->rtol = rtol; 1199 } 1200 1201 if (atol != (PetscReal)PETSC_DEFAULT) { 1202 PetscCheck(atol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Absolute tolerance %14.12e must be non-negative", (double)atol); 1203 linesearch->atol = atol; 1204 } 1205 1206 if (ltol != (PetscReal)PETSC_DEFAULT) { 1207 PetscCheck(ltol >= 0.0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Lambda tolerance %14.12e must be non-negative", (double)ltol); 1208 linesearch->ltol = ltol; 1209 } 1210 1211 if (max_it != PETSC_DEFAULT) { 1212 PetscCheck(max_it >= 0, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_ARG_OUTOFRANGE, "Maximum number of iterations %" PetscInt_FMT " must be non-negative", max_it); 1213 linesearch->max_its = max_it; 1214 } 1215 PetscFunctionReturn(PETSC_SUCCESS); 1216 } 1217 1218 /*@ 1219 SNESLineSearchGetDamping - Gets the line search damping parameter. 1220 1221 Input Parameter: 1222 . linesearch - the line search context 1223 1224 Output Parameter: 1225 . damping - The damping parameter 1226 1227 Level: advanced 1228 1229 .seealso: [](ch_snes), `SNES`, `SNESLineSearchGetStepTolerance()`, `SNESQN` 1230 @*/ 1231 PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch, PetscReal *damping) 1232 { 1233 PetscFunctionBegin; 1234 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1235 PetscAssertPointer(damping, 2); 1236 *damping = linesearch->damping; 1237 PetscFunctionReturn(PETSC_SUCCESS); 1238 } 1239 1240 /*@ 1241 SNESLineSearchSetDamping - Sets the line search damping parameter. 1242 1243 Input Parameters: 1244 + linesearch - the line search context 1245 - damping - The damping parameter 1246 1247 Options Database Key: 1248 . -snes_linesearch_damping <damping> - the damping value 1249 1250 Level: intermediate 1251 1252 Note: 1253 The `SNESLINESEARCHNONE` line search merely takes the update step scaled by the damping parameter. 1254 The use of the damping parameter in the `SNESLINESEARCHL2` and `SNESLINESEARCHCP` line searches is much more subtle; 1255 it is used as a starting point in calculating the secant step. However, the eventual 1256 step may be of greater length than the damping parameter. In the `SNESLINESEARCHBT` line search it is 1257 used as the maximum possible step length, as the `SNESLINESEARCHBT` line search only backtracks. 1258 1259 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetDamping()` 1260 @*/ 1261 PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch, PetscReal damping) 1262 { 1263 PetscFunctionBegin; 1264 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1265 linesearch->damping = damping; 1266 PetscFunctionReturn(PETSC_SUCCESS); 1267 } 1268 1269 /*@ 1270 SNESLineSearchGetOrder - Gets the line search approximation order. 1271 1272 Input Parameter: 1273 . linesearch - the line search context 1274 1275 Output Parameter: 1276 . order - The order 1277 1278 Level: intermediate 1279 1280 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetOrder()` 1281 @*/ 1282 PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch, PetscInt *order) 1283 { 1284 PetscFunctionBegin; 1285 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1286 PetscAssertPointer(order, 2); 1287 *order = linesearch->order; 1288 PetscFunctionReturn(PETSC_SUCCESS); 1289 } 1290 1291 /*@ 1292 SNESLineSearchSetOrder - Sets the maximum order of the polynomial fit used in the line search 1293 1294 Input Parameters: 1295 + linesearch - the line search context 1296 - order - The order 1297 1298 Level: intermediate 1299 1300 Values for `order`\: 1301 + 1 or `SNES_LINESEARCH_ORDER_LINEAR` - linear order 1302 . 2 or `SNES_LINESEARCH_ORDER_QUADRATIC` - quadratic order 1303 - 3 or `SNES_LINESEARCH_ORDER_CUBIC` - cubic order 1304 1305 Options Database Key: 1306 . -snes_linesearch_order <order> - 1, 2, 3. Most types only support certain orders (`SNESLINESEARCHBT` supports 2 or 3) 1307 1308 Note: 1309 These orders are supported by `SNESLINESEARCHBT` and `SNESLINESEARCHCP` 1310 1311 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetOrder()`, `SNESLineSearchSetDamping()` 1312 @*/ 1313 PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch, PetscInt order) 1314 { 1315 PetscFunctionBegin; 1316 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1317 linesearch->order = order; 1318 PetscFunctionReturn(PETSC_SUCCESS); 1319 } 1320 1321 /*@ 1322 SNESLineSearchGetNorms - Gets the norms for the current solution `X`, the current update `Y`, and the current function value `F`. 1323 1324 Not Collective 1325 1326 Input Parameter: 1327 . linesearch - the line search context 1328 1329 Output Parameters: 1330 + xnorm - The norm of the current solution 1331 . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution. 1332 - ynorm - The norm of the current update (after scaling by the linesearch computed lambda) 1333 1334 Level: developer 1335 1336 Notes: 1337 Some values may not be up-to-date at particular points in the code. 1338 1339 This, in combination with `SNESLineSearchSetNorms()`, allow the line search and the `SNESSolve_XXX()` to share 1340 computed values. 1341 1342 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()` 1343 @*/ 1344 PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal *xnorm, PetscReal *fnorm, PetscReal *ynorm) 1345 { 1346 PetscFunctionBegin; 1347 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1348 if (xnorm) *xnorm = linesearch->xnorm; 1349 if (fnorm) *fnorm = linesearch->fnorm; 1350 if (ynorm) *ynorm = linesearch->ynorm; 1351 PetscFunctionReturn(PETSC_SUCCESS); 1352 } 1353 1354 /*@ 1355 SNESLineSearchSetNorms - Sets the computed norms for the current solution `X`, the current update `Y`, and the current function value `F`. 1356 1357 Collective 1358 1359 Input Parameters: 1360 + linesearch - the line search context 1361 . xnorm - The norm of the current solution 1362 . fnorm - The norm of the current function, this is the `norm(function(X))` where `X` is the current solution 1363 - ynorm - The norm of the current update (after scaling by the linesearch computed lambda) 1364 1365 Level: developer 1366 1367 Note: 1368 This is called by the line search routines to store the values they have just computed 1369 1370 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1371 @*/ 1372 PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 1373 { 1374 PetscFunctionBegin; 1375 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1376 linesearch->xnorm = xnorm; 1377 linesearch->fnorm = fnorm; 1378 linesearch->ynorm = ynorm; 1379 PetscFunctionReturn(PETSC_SUCCESS); 1380 } 1381 1382 /*@ 1383 SNESLineSearchComputeNorms - Explicitly computes the norms of the current solution `X`, the current update `Y`, and the current function value `F`. 1384 1385 Input Parameter: 1386 . linesearch - the line search context 1387 1388 Options Database Key: 1389 . -snes_linesearch_norms - turn norm computation on or off 1390 1391 Level: intermediate 1392 1393 Developer Note: 1394 The options database key is misnamed. It should be -snes_linesearch_compute_norms 1395 1396 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms`, `SNESLineSearchSetNorms()`, `SNESLineSearchSetComputeNorms()` 1397 @*/ 1398 PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 1399 { 1400 SNES snes; 1401 1402 PetscFunctionBegin; 1403 if (linesearch->norms) { 1404 if (linesearch->ops->vinorm) { 1405 PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 1406 PetscCall(VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 1407 PetscCall(VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 1408 PetscCall((*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm)); 1409 } else { 1410 PetscCall(VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 1411 PetscCall(VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 1412 PetscCall(VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 1413 PetscCall(VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm)); 1414 PetscCall(VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm)); 1415 PetscCall(VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm)); 1416 } 1417 } 1418 PetscFunctionReturn(PETSC_SUCCESS); 1419 } 1420 1421 /*@ 1422 SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 1423 1424 Input Parameters: 1425 + linesearch - the line search context 1426 - flg - indicates whether or not to compute norms 1427 1428 Options Database Key: 1429 . -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic (none) `SNESLINESEARCHBASIC` line search 1430 1431 Level: intermediate 1432 1433 Note: 1434 This is most relevant to the `SNESLINESEARCHBASIC` (or equivalently `SNESLINESEARCHNONE`) line search type since most line searches have a stopping criteria involving the norm. 1435 1436 Developer Note: 1437 The options database key is misnamed. It should be -snes_linesearch_compute_norms 1438 1439 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetNorms()`, `SNESLineSearchComputeNorms()`, `SNESLINESEARCHBASIC` 1440 @*/ 1441 PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 1442 { 1443 PetscFunctionBegin; 1444 linesearch->norms = flg; 1445 PetscFunctionReturn(PETSC_SUCCESS); 1446 } 1447 1448 /*@ 1449 SNESLineSearchGetVecs - Gets the vectors from the `SNESLineSearch` context 1450 1451 Not Collective but the vectors are parallel 1452 1453 Input Parameter: 1454 . linesearch - the line search context 1455 1456 Output Parameters: 1457 + X - Solution vector 1458 . F - Function vector 1459 . Y - Search direction vector 1460 . W - Solution work vector 1461 - G - Function work vector 1462 1463 Level: advanced 1464 1465 Notes: 1466 At the beginning of a line search application, `X` should contain a 1467 solution and the vector `F` the function computed at `X`. At the end of the 1468 line search application, `X` should contain the new solution, and `F` the 1469 function evaluated at the new solution. 1470 1471 These vectors are owned by the `SNESLineSearch` and should not be destroyed by the caller 1472 1473 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetNorms()`, `SNESLineSearchSetVecs()` 1474 @*/ 1475 PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch, Vec *X, Vec *F, Vec *Y, Vec *W, Vec *G) 1476 { 1477 PetscFunctionBegin; 1478 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1479 if (X) { 1480 PetscAssertPointer(X, 2); 1481 *X = linesearch->vec_sol; 1482 } 1483 if (F) { 1484 PetscAssertPointer(F, 3); 1485 *F = linesearch->vec_func; 1486 } 1487 if (Y) { 1488 PetscAssertPointer(Y, 4); 1489 *Y = linesearch->vec_update; 1490 } 1491 if (W) { 1492 PetscAssertPointer(W, 5); 1493 *W = linesearch->vec_sol_new; 1494 } 1495 if (G) { 1496 PetscAssertPointer(G, 6); 1497 *G = linesearch->vec_func_new; 1498 } 1499 PetscFunctionReturn(PETSC_SUCCESS); 1500 } 1501 1502 /*@ 1503 SNESLineSearchSetVecs - Sets the vectors on the `SNESLineSearch` context 1504 1505 Logically Collective 1506 1507 Input Parameters: 1508 + linesearch - the line search context 1509 . X - Solution vector 1510 . F - Function vector 1511 . Y - Search direction vector 1512 . W - Solution work vector 1513 - G - Function work vector 1514 1515 Level: developer 1516 1517 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetNorms()`, `SNESLineSearchGetVecs()` 1518 @*/ 1519 PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch, Vec X, Vec F, Vec Y, Vec W, Vec G) 1520 { 1521 PetscFunctionBegin; 1522 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1523 if (X) { 1524 PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 1525 linesearch->vec_sol = X; 1526 } 1527 if (F) { 1528 PetscValidHeaderSpecific(F, VEC_CLASSID, 3); 1529 linesearch->vec_func = F; 1530 } 1531 if (Y) { 1532 PetscValidHeaderSpecific(Y, VEC_CLASSID, 4); 1533 linesearch->vec_update = Y; 1534 } 1535 if (W) { 1536 PetscValidHeaderSpecific(W, VEC_CLASSID, 5); 1537 linesearch->vec_sol_new = W; 1538 } 1539 if (G) { 1540 PetscValidHeaderSpecific(G, VEC_CLASSID, 6); 1541 linesearch->vec_func_new = G; 1542 } 1543 PetscFunctionReturn(PETSC_SUCCESS); 1544 } 1545 1546 /*@ 1547 SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1548 `SNESLineSearch` options in the database. 1549 1550 Logically Collective 1551 1552 Input Parameters: 1553 + linesearch - the `SNESLineSearch` context 1554 - prefix - the prefix to prepend to all option names 1555 1556 Level: advanced 1557 1558 Note: 1559 A hyphen (-) must NOT be given at the beginning of the prefix name. 1560 The first character of all runtime options is AUTOMATICALLY the hyphen. 1561 1562 .seealso: [](ch_snes), `SNES`, `SNESLineSearch()`, `SNESLineSearchSetFromOptions()`, `SNESGetOptionsPrefix()` 1563 @*/ 1564 PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch, const char prefix[]) 1565 { 1566 PetscFunctionBegin; 1567 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1568 PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)linesearch, prefix)); 1569 PetscFunctionReturn(PETSC_SUCCESS); 1570 } 1571 1572 /*@ 1573 SNESLineSearchGetOptionsPrefix - Gets the prefix used for searching for all 1574 SNESLineSearch options in the database. 1575 1576 Not Collective 1577 1578 Input Parameter: 1579 . linesearch - the `SNESLineSearch` context 1580 1581 Output Parameter: 1582 . prefix - pointer to the prefix string used 1583 1584 Level: advanced 1585 1586 Fortran Notes: 1587 The user should pass in a string 'prefix' of 1588 sufficient length to hold the prefix. 1589 1590 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESAppendOptionsPrefix()` 1591 @*/ 1592 PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch, const char *prefix[]) 1593 { 1594 PetscFunctionBegin; 1595 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1596 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)linesearch, prefix)); 1597 PetscFunctionReturn(PETSC_SUCCESS); 1598 } 1599 1600 /*@C 1601 SNESLineSearchSetWorkVecs - Sets work vectors for the line search. 1602 1603 Input Parameters: 1604 + linesearch - the `SNESLineSearch` context 1605 - nwork - the number of work vectors 1606 1607 Level: developer 1608 1609 Developer Note: 1610 This is called from within the set up routines for each of the line search types `SNESLineSearchType` 1611 1612 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESSetWorkVecs()` 1613 @*/ 1614 PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) 1615 { 1616 PetscFunctionBegin; 1617 PetscCheck(linesearch->vec_sol, PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 1618 PetscCall(VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work)); 1619 PetscFunctionReturn(PETSC_SUCCESS); 1620 } 1621 1622 /*@ 1623 SNESLineSearchGetReason - Gets the success/failure status of the last line search application 1624 1625 Input Parameter: 1626 . linesearch - the line search context 1627 1628 Output Parameter: 1629 . result - The success or failure status 1630 1631 Level: developer 1632 1633 Note: 1634 This is typically called after `SNESLineSearchApply()` in order to determine if the line search failed 1635 (and set into the `SNES` convergence accordingly). 1636 1637 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetReason()`, `SNESLineSearchReason` 1638 @*/ 1639 PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) 1640 { 1641 PetscFunctionBegin; 1642 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1643 PetscAssertPointer(result, 2); 1644 *result = linesearch->result; 1645 PetscFunctionReturn(PETSC_SUCCESS); 1646 } 1647 1648 /*@ 1649 SNESLineSearchSetReason - Sets the success/failure status of the line search application 1650 1651 Logically Collective; No Fortran Support 1652 1653 Input Parameters: 1654 + linesearch - the line search context 1655 - result - The success or failure status 1656 1657 Level: developer 1658 1659 Note: 1660 This is typically called in a `SNESLineSearchType` implementation of `SNESLineSearchApply()` or a `SNESLINESEARCHSHELL` implementation to set 1661 the success or failure of the line search method. 1662 1663 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchReason`, `SNESLineSearchGetSResult()` 1664 @*/ 1665 PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) 1666 { 1667 PetscFunctionBegin; 1668 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1669 linesearch->result = result; 1670 PetscFunctionReturn(PETSC_SUCCESS); 1671 } 1672 1673 // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation 1674 /*@C 1675 SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 1676 1677 Logically Collective 1678 1679 Input Parameters: 1680 + linesearch - the linesearch object 1681 . projectfunc - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence 1682 - normfunc - function for computing the norm of an active set, see `SNESLineSearchVINormFn` for calling sequence 1683 1684 Level: advanced 1685 1686 Notes: 1687 The VI solvers require projection of the solution to the feasible set. `projectfunc` should implement this. 1688 1689 The VI solvers require special evaluation of the function norm such that the norm is only calculated 1690 on the inactive set. This should be implemented by `normfunc`. 1691 1692 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchGetVIFunctions()`, `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()`, 1693 `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn` 1694 @*/ 1695 PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn *projectfunc, SNESLineSearchVINormFn *normfunc) 1696 { 1697 PetscFunctionBegin; 1698 PetscValidHeaderSpecific(linesearch, SNESLINESEARCH_CLASSID, 1); 1699 if (projectfunc) linesearch->ops->viproject = projectfunc; 1700 if (normfunc) linesearch->ops->vinorm = normfunc; 1701 PetscFunctionReturn(PETSC_SUCCESS); 1702 } 1703 1704 /*@C 1705 SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 1706 1707 Not Collective 1708 1709 Input Parameter: 1710 . linesearch - the line search context, obtain with `SNESGetLineSearch()` 1711 1712 Output Parameters: 1713 + projectfunc - function for projecting the function to the bounds, see `SNESLineSearchVIProjectFn` for calling sequence 1714 - normfunc - function for computing the norm of an active set, see `SNESLineSearchVINormFn ` for calling sequence 1715 1716 Level: advanced 1717 1718 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchSetVIFunctions()`, `SNESLineSearchGetPostCheck()`, `SNESLineSearchGetPreCheck()`, 1719 `SNESLineSearchVIProjectFn`, `SNESLineSearchVINormFn` 1720 @*/ 1721 PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFn **projectfunc, SNESLineSearchVINormFn **normfunc) 1722 { 1723 PetscFunctionBegin; 1724 if (projectfunc) *projectfunc = linesearch->ops->viproject; 1725 if (normfunc) *normfunc = linesearch->ops->vinorm; 1726 PetscFunctionReturn(PETSC_SUCCESS); 1727 } 1728 1729 /*@C 1730 SNESLineSearchRegister - register a line search type `SNESLineSearchType` 1731 1732 Logically Collective, No Fortran Support 1733 1734 Input Parameters: 1735 + sname - name of the `SNESLineSearchType()` 1736 - function - the creation function for that type 1737 1738 Calling sequence of `function`: 1739 . ls - the line search context 1740 1741 Level: advanced 1742 1743 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchSetType()` 1744 @*/ 1745 PetscErrorCode SNESLineSearchRegister(const char sname[], PetscErrorCode (*function)(SNESLineSearch ls)) 1746 { 1747 PetscFunctionBegin; 1748 PetscCall(SNESInitializePackage()); 1749 PetscCall(PetscFunctionListAdd(&SNESLineSearchList, sname, function)); 1750 PetscFunctionReturn(PETSC_SUCCESS); 1751 } 1752