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