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