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