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