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