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(PetscOptionsObject,(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 Options Database: 1133 . -snes_linesearch_damping 1134 Level: intermediate 1135 1136 Notes: 1137 The basic line search merely takes the update step scaled by the damping parameter. 1138 The use of the damping parameter in the l2 and cp line searches is much more subtle; 1139 it is used as a starting point in calculating the secant step. However, the eventual 1140 step may be of greater length than the damping parameter. In the bt line search it is 1141 used as the maximum possible step length, as the bt line search only backtracks. 1142 1143 .seealso: SNESLineSearchGetDamping() 1144 @*/ 1145 PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping) 1146 { 1147 PetscFunctionBegin; 1148 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1149 linesearch->damping = damping; 1150 PetscFunctionReturn(0); 1151 } 1152 1153 #undef __FUNCT__ 1154 #define __FUNCT__ "SNESLineSearchGetOrder" 1155 /*@ 1156 SNESLineSearchGetOrder - Gets the line search approximation order. 1157 1158 Input Parameters: 1159 . linesearch - linesearch context 1160 1161 Output Parameters: 1162 . order - The order 1163 1164 Possible Values for order: 1165 + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 1166 . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 1167 - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 1168 1169 Level: intermediate 1170 1171 .seealso: SNESLineSearchSetOrder() 1172 @*/ 1173 1174 PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order) 1175 { 1176 PetscFunctionBegin; 1177 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1178 PetscValidPointer(order, 2); 1179 *order = linesearch->order; 1180 PetscFunctionReturn(0); 1181 } 1182 1183 #undef __FUNCT__ 1184 #define __FUNCT__ "SNESLineSearchSetOrder" 1185 /*@ 1186 SNESLineSearchSetOrder - Sets the line search damping paramter. 1187 1188 Input Parameters: 1189 . linesearch - linesearch context 1190 . order - The damping parameter 1191 1192 Level: intermediate 1193 1194 Possible Values for order: 1195 + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 1196 . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 1197 - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 1198 1199 Notes: 1200 Variable orders are supported by the following line searches: 1201 + bt - cubic and quadratic 1202 - cp - linear and quadratic 1203 1204 .seealso: SNESLineSearchGetOrder() 1205 @*/ 1206 PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order) 1207 { 1208 PetscFunctionBegin; 1209 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1210 linesearch->order = order; 1211 PetscFunctionReturn(0); 1212 } 1213 1214 #undef __FUNCT__ 1215 #define __FUNCT__ "SNESLineSearchGetNorms" 1216 /*@ 1217 SNESLineSearchGetNorms - Gets the norms for for X, Y, and F. 1218 1219 Input Parameters: 1220 . linesearch - linesearch context 1221 1222 Output Parameters: 1223 + xnorm - The norm of the current solution 1224 . fnorm - The norm of the current function 1225 - ynorm - The norm of the current update 1226 1227 Notes: 1228 This function is mainly called from SNES implementations. 1229 1230 Level: developer 1231 1232 .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs() 1233 @*/ 1234 PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm) 1235 { 1236 PetscFunctionBegin; 1237 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1238 if (xnorm) *xnorm = linesearch->xnorm; 1239 if (fnorm) *fnorm = linesearch->fnorm; 1240 if (ynorm) *ynorm = linesearch->ynorm; 1241 PetscFunctionReturn(0); 1242 } 1243 1244 #undef __FUNCT__ 1245 #define __FUNCT__ "SNESLineSearchSetNorms" 1246 /*@ 1247 SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 1248 1249 Input Parameters: 1250 + linesearch - linesearch context 1251 . xnorm - The norm of the current solution 1252 . fnorm - The norm of the current function 1253 - ynorm - The norm of the current update 1254 1255 Level: advanced 1256 1257 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1258 @*/ 1259 PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 1260 { 1261 PetscFunctionBegin; 1262 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1263 linesearch->xnorm = xnorm; 1264 linesearch->fnorm = fnorm; 1265 linesearch->ynorm = ynorm; 1266 PetscFunctionReturn(0); 1267 } 1268 1269 #undef __FUNCT__ 1270 #define __FUNCT__ "SNESLineSearchComputeNorms" 1271 /*@ 1272 SNESLineSearchComputeNorms - Computes the norms of X, F, and Y. 1273 1274 Input Parameters: 1275 . linesearch - linesearch context 1276 1277 Options Database Keys: 1278 . -snes_linesearch_norms - turn norm computation on or off 1279 1280 Level: intermediate 1281 1282 .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms() 1283 @*/ 1284 PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 1285 { 1286 PetscErrorCode ierr; 1287 SNES snes; 1288 1289 PetscFunctionBegin; 1290 if (linesearch->norms) { 1291 if (linesearch->ops->vinorm) { 1292 ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 1293 ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1294 ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1295 ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr); 1296 } else { 1297 ierr = VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 1298 ierr = VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1299 ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1300 ierr = VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 1301 ierr = VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1302 ierr = VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1303 } 1304 } 1305 PetscFunctionReturn(0); 1306 } 1307 1308 #undef __FUNCT__ 1309 #define __FUNCT__ "SNESLineSearchSetComputeNorms" 1310 /*@ 1311 SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 1312 1313 Input Parameters: 1314 + linesearch - linesearch context 1315 - flg - indicates whether or not to compute norms 1316 1317 Options Database Keys: 1318 . -snes_linesearch_norms - turn norm computation on or off 1319 1320 Notes: 1321 This is most relevant to the SNESLINESEARCHBASIC line search type. 1322 1323 Level: intermediate 1324 1325 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC 1326 @*/ 1327 PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 1328 { 1329 PetscFunctionBegin; 1330 linesearch->norms = flg; 1331 PetscFunctionReturn(0); 1332 } 1333 1334 #undef __FUNCT__ 1335 #define __FUNCT__ "SNESLineSearchGetVecs" 1336 /*@ 1337 SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context 1338 1339 Input Parameters: 1340 . linesearch - linesearch context 1341 1342 Output Parameters: 1343 + X - Solution vector 1344 . F - Function vector 1345 . Y - Search direction vector 1346 . W - Solution work vector 1347 - G - Function work vector 1348 1349 Notes: 1350 At the beginning of a line search application, X should contain a 1351 solution and the vector F the function computed at X. At the end of the 1352 line search application, X should contain the new solution, and F the 1353 function evaluated at the new solution. 1354 1355 These vectors are owned by the SNESLineSearch and should not be destroyed by the caller 1356 1357 Level: advanced 1358 1359 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1360 @*/ 1361 PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) 1362 { 1363 PetscFunctionBegin; 1364 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1365 if (X) { 1366 PetscValidPointer(X, 2); 1367 *X = linesearch->vec_sol; 1368 } 1369 if (F) { 1370 PetscValidPointer(F, 3); 1371 *F = linesearch->vec_func; 1372 } 1373 if (Y) { 1374 PetscValidPointer(Y, 4); 1375 *Y = linesearch->vec_update; 1376 } 1377 if (W) { 1378 PetscValidPointer(W, 5); 1379 *W = linesearch->vec_sol_new; 1380 } 1381 if (G) { 1382 PetscValidPointer(G, 6); 1383 *G = linesearch->vec_func_new; 1384 } 1385 PetscFunctionReturn(0); 1386 } 1387 1388 #undef __FUNCT__ 1389 #define __FUNCT__ "SNESLineSearchSetVecs" 1390 /*@ 1391 SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context 1392 1393 Input Parameters: 1394 + linesearch - linesearch context 1395 . X - Solution vector 1396 . F - Function vector 1397 . Y - Search direction vector 1398 . W - Solution work vector 1399 - G - Function work vector 1400 1401 Level: advanced 1402 1403 .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs() 1404 @*/ 1405 PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) 1406 { 1407 PetscFunctionBegin; 1408 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1409 if (X) { 1410 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1411 linesearch->vec_sol = X; 1412 } 1413 if (F) { 1414 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 1415 linesearch->vec_func = F; 1416 } 1417 if (Y) { 1418 PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 1419 linesearch->vec_update = Y; 1420 } 1421 if (W) { 1422 PetscValidHeaderSpecific(W,VEC_CLASSID,5); 1423 linesearch->vec_sol_new = W; 1424 } 1425 if (G) { 1426 PetscValidHeaderSpecific(G,VEC_CLASSID,6); 1427 linesearch->vec_func_new = G; 1428 } 1429 PetscFunctionReturn(0); 1430 } 1431 1432 #undef __FUNCT__ 1433 #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix" 1434 /*@C 1435 SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1436 SNES options in the database. 1437 1438 Logically Collective on SNESLineSearch 1439 1440 Input Parameters: 1441 + snes - the SNES context 1442 - prefix - the prefix to prepend to all option names 1443 1444 Notes: 1445 A hyphen (-) must NOT be given at the beginning of the prefix name. 1446 The first character of all runtime options is AUTOMATICALLY the hyphen. 1447 1448 Level: advanced 1449 1450 .keywords: SNESLineSearch, append, options, prefix, database 1451 1452 .seealso: SNESGetOptionsPrefix() 1453 @*/ 1454 PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[]) 1455 { 1456 PetscErrorCode ierr; 1457 1458 PetscFunctionBegin; 1459 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1460 ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1461 PetscFunctionReturn(0); 1462 } 1463 1464 #undef __FUNCT__ 1465 #define __FUNCT__ "SNESLineSearchGetOptionsPrefix" 1466 /*@C 1467 SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all 1468 SNESLineSearch options in the database. 1469 1470 Not Collective 1471 1472 Input Parameter: 1473 . linesearch - the SNESLineSearch context 1474 1475 Output Parameter: 1476 . prefix - pointer to the prefix string used 1477 1478 Notes: 1479 On the fortran side, the user should pass in a string 'prefix' of 1480 sufficient length to hold the prefix. 1481 1482 Level: advanced 1483 1484 .keywords: SNESLineSearch, get, options, prefix, database 1485 1486 .seealso: SNESAppendOptionsPrefix() 1487 @*/ 1488 PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[]) 1489 { 1490 PetscErrorCode ierr; 1491 1492 PetscFunctionBegin; 1493 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1494 ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1495 PetscFunctionReturn(0); 1496 } 1497 1498 #undef __FUNCT__ 1499 #define __FUNCT__ "SNESLineSearchSetWorkVecs" 1500 /*@C 1501 SNESLineSearchSetWorkVecs - Gets work vectors for the line search. 1502 1503 Input Parameter: 1504 + linesearch - the SNESLineSearch context 1505 - nwork - the number of work vectors 1506 1507 Level: developer 1508 1509 Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations 1510 1511 .keywords: SNESLineSearch, work, vector 1512 1513 .seealso: SNESSetWorkVecs() 1514 @*/ 1515 PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) 1516 { 1517 PetscErrorCode ierr; 1518 1519 PetscFunctionBegin; 1520 if (linesearch->vec_sol) { 1521 ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr); 1522 } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 1523 PetscFunctionReturn(0); 1524 } 1525 1526 #undef __FUNCT__ 1527 #define __FUNCT__ "SNESLineSearchGetReason" 1528 /*@ 1529 SNESLineSearchGetReason - Gets the success/failure status of the last line search application 1530 1531 Input Parameters: 1532 . linesearch - linesearch context 1533 1534 Output Parameters: 1535 . result - The success or failure status 1536 1537 Notes: 1538 This is typically called after SNESLineSearchApply() in order to determine if the line-search failed 1539 (and set the SNES convergence accordingly). 1540 1541 Level: intermediate 1542 1543 .seealso: SNESLineSearchSetReason(), SNESLineSearchReason 1544 @*/ 1545 PetscErrorCode SNESLineSearchGetReason(SNESLineSearch linesearch, SNESLineSearchReason *result) 1546 { 1547 PetscFunctionBegin; 1548 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1549 PetscValidPointer(result, 2); 1550 *result = linesearch->result; 1551 PetscFunctionReturn(0); 1552 } 1553 1554 #undef __FUNCT__ 1555 #define __FUNCT__ "SNESLineSearchSetReason" 1556 /*@ 1557 SNESLineSearchSetReason - Sets the success/failure status of the last line search application 1558 1559 Input Parameters: 1560 + linesearch - linesearch context 1561 - result - The success or failure status 1562 1563 Notes: 1564 This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set 1565 the success or failure of the line search method. 1566 1567 Level: developer 1568 1569 .seealso: SNESLineSearchGetSResult() 1570 @*/ 1571 PetscErrorCode SNESLineSearchSetReason(SNESLineSearch linesearch, SNESLineSearchReason result) 1572 { 1573 PetscFunctionBegin; 1574 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1575 linesearch->result = result; 1576 PetscFunctionReturn(0); 1577 } 1578 1579 #undef __FUNCT__ 1580 #define __FUNCT__ "SNESLineSearchSetVIFunctions" 1581 /*@C 1582 SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 1583 1584 Input Parameters: 1585 + snes - nonlinear context obtained from SNESCreate() 1586 . projectfunc - function for projecting the function to the bounds 1587 - normfunc - function for computing the norm of an active set 1588 1589 Logically Collective on SNES 1590 1591 Calling sequence of projectfunc: 1592 .vb 1593 projectfunc (SNES snes, Vec X) 1594 .ve 1595 1596 Input parameters for projectfunc: 1597 + snes - nonlinear context 1598 - X - current solution 1599 1600 Output parameters for projectfunc: 1601 . X - Projected solution 1602 1603 Calling sequence of normfunc: 1604 .vb 1605 projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm) 1606 .ve 1607 1608 Input parameters for normfunc: 1609 + snes - nonlinear context 1610 . X - current solution 1611 - F - current residual 1612 1613 Output parameters for normfunc: 1614 . fnorm - VI-specific norm of the function 1615 1616 Notes: 1617 The VI solvers require projection of the solution to the feasible set. projectfunc should implement this. 1618 1619 The VI solvers require special evaluation of the function norm such that the norm is only calculated 1620 on the inactive set. This should be implemented by normfunc. 1621 1622 Level: developer 1623 1624 .keywords: SNES, line search, VI, nonlinear, set, line search 1625 1626 .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() 1627 @*/ 1628 extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 1629 { 1630 PetscFunctionBegin; 1631 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1632 if (projectfunc) linesearch->ops->viproject = projectfunc; 1633 if (normfunc) linesearch->ops->vinorm = normfunc; 1634 PetscFunctionReturn(0); 1635 } 1636 1637 #undef __FUNCT__ 1638 #define __FUNCT__ "SNESLineSearchGetVIFunctions" 1639 /*@C 1640 SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 1641 1642 Input Parameters: 1643 . linesearch - the line search context, obtain with SNESGetLineSearch() 1644 1645 Output Parameters: 1646 + projectfunc - function for projecting the function to the bounds 1647 - normfunc - function for computing the norm of an active set 1648 1649 Logically Collective on SNES 1650 1651 Level: developer 1652 1653 .keywords: SNES, line search, VI, nonlinear, get, line search 1654 1655 .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck() 1656 @*/ 1657 extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 1658 { 1659 PetscFunctionBegin; 1660 if (projectfunc) *projectfunc = linesearch->ops->viproject; 1661 if (normfunc) *normfunc = linesearch->ops->vinorm; 1662 PetscFunctionReturn(0); 1663 } 1664 1665 #undef __FUNCT__ 1666 #define __FUNCT__ "SNESLineSearchRegister" 1667 /*@C 1668 SNESLineSearchRegister - See SNESLineSearchRegister() 1669 1670 Level: advanced 1671 @*/ 1672 PetscErrorCode SNESLineSearchRegister(const char sname[],PetscErrorCode (*function)(SNESLineSearch)) 1673 { 1674 PetscErrorCode ierr; 1675 1676 PetscFunctionBegin; 1677 ierr = PetscFunctionListAdd(&SNESLineSearchList,sname,function);CHKERRQ(ierr); 1678 PetscFunctionReturn(0); 1679 } 1680