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