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 = PetscObjectAMSUnPublish((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 if (PetscAMSPublishAll) { 823 ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr); 824 } 825 PetscFunctionReturn(0); 826 } 827 828 #undef __FUNCT__ 829 #define __FUNCT__ "SNESLineSearchSetSNES" 830 /*@ 831 SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation. 832 833 Input Parameters: 834 + linesearch - linesearch context 835 - snes - The snes instance 836 837 Level: developer 838 839 Notes: 840 This happens automatically when the line search is gotten/created with 841 SNESGetSNESLineSearch(). This routine is therefore mainly called within SNES 842 implementations. 843 844 Level: developer 845 846 .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 847 @*/ 848 PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes) 849 { 850 PetscFunctionBegin; 851 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 852 PetscValidHeaderSpecific(snes,SNES_CLASSID,2); 853 linesearch->snes = snes; 854 PetscFunctionReturn(0); 855 } 856 857 #undef __FUNCT__ 858 #define __FUNCT__ "SNESLineSearchGetSNES" 859 /*@ 860 SNESLineSearchGetSNES - Gets the SNES instance associated with the line search. 861 Having an associated SNES is necessary because most line search implementations must be able to 862 evaluate the function using SNESComputeFunction() for the associated SNES. This routine 863 is used in the line search implementations when one must get this associated SNES instance. 864 865 Input Parameters: 866 . linesearch - linesearch context 867 868 Output Parameters: 869 . snes - The snes instance 870 871 Level: developer 872 873 .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 874 @*/ 875 PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes) 876 { 877 PetscFunctionBegin; 878 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 879 PetscValidPointer(snes, 2); 880 *snes = linesearch->snes; 881 PetscFunctionReturn(0); 882 } 883 884 #undef __FUNCT__ 885 #define __FUNCT__ "SNESLineSearchGetLambda" 886 /*@ 887 SNESLineSearchGetLambda - Gets the last linesearch steplength discovered. 888 889 Input Parameters: 890 . linesearch - linesearch context 891 892 Output Parameters: 893 . lambda - The last steplength computed during SNESLineSearchApply() 894 895 Level: advanced 896 897 Notes: 898 This is useful in methods where the solver is ill-scaled and 899 requires some adaptive notion of the difference in scale between the 900 solution and the function. For instance, SNESQN may be scaled by the 901 line search lambda using the argument -snes_qn_scaling ls. 902 903 904 .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply() 905 @*/ 906 PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda) 907 { 908 PetscFunctionBegin; 909 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 910 PetscValidPointer(lambda, 2); 911 *lambda = linesearch->lambda; 912 PetscFunctionReturn(0); 913 } 914 915 #undef __FUNCT__ 916 #define __FUNCT__ "SNESLineSearchSetLambda" 917 /*@ 918 SNESLineSearchSetLambda - Sets the linesearch steplength. 919 920 Input Parameters: 921 + linesearch - linesearch context 922 - lambda - The last steplength. 923 924 Notes: 925 This routine is typically used within implementations of SNESLineSearchApply 926 to set the final steplength. This routine (and SNESLineSearchGetLambda()) were 927 added in order to facilitate Quasi-Newton methods that use the previous steplength 928 as an inner scaling parameter. 929 930 Level: advanced 931 932 .seealso: SNESLineSearchGetLambda() 933 @*/ 934 PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 935 { 936 PetscFunctionBegin; 937 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 938 linesearch->lambda = lambda; 939 PetscFunctionReturn(0); 940 } 941 942 #undef __FUNCT__ 943 #define __FUNCT__ "SNESLineSearchGetTolerances" 944 /*@ 945 SNESLineSearchGetTolerances - Gets the tolerances for the linesearch. These include 946 tolerances for the relative and absolute change in the function norm, the change 947 in lambda for iterative line searches, the minimum steplength, the maximum steplength, 948 and the maximum number of iterations the line search procedure may take. 949 950 Input Parameters: 951 . linesearch - linesearch context 952 953 Output Parameters: 954 + steptol - The minimum steplength 955 . maxstep - The maximum steplength 956 . rtol - The relative tolerance for iterative line searches 957 . atol - The absolute tolerance for iterative line searches 958 . ltol - The change in lambda tolerance for iterative line searches 959 - max_it - The maximum number of iterations of the line search 960 961 Level: intermediate 962 963 Notes: 964 Different line searches may implement these parameters slightly differently as 965 the type requires. 966 967 .seealso: SNESLineSearchSetTolerances() 968 @*/ 969 PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 970 { 971 PetscFunctionBegin; 972 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 973 if (steptol) { 974 PetscValidPointer(steptol, 2); 975 *steptol = linesearch->steptol; 976 } 977 if (maxstep) { 978 PetscValidPointer(maxstep, 3); 979 *maxstep = linesearch->maxstep; 980 } 981 if (rtol) { 982 PetscValidPointer(rtol, 4); 983 *rtol = linesearch->rtol; 984 } 985 if (atol) { 986 PetscValidPointer(atol, 5); 987 *atol = linesearch->atol; 988 } 989 if (ltol) { 990 PetscValidPointer(ltol, 6); 991 *ltol = linesearch->ltol; 992 } 993 if (max_its) { 994 PetscValidPointer(max_its, 7); 995 *max_its = linesearch->max_its; 996 } 997 PetscFunctionReturn(0); 998 } 999 1000 #undef __FUNCT__ 1001 #define __FUNCT__ "SNESLineSearchSetTolerances" 1002 /*@ 1003 SNESLineSearchSetTolerances - Gets the tolerances for the linesearch. These include 1004 tolerances for the relative and absolute change in the function norm, the change 1005 in lambda for iterative line searches, the minimum steplength, the maximum steplength, 1006 and the maximum number of iterations the line search procedure may take. 1007 1008 Input Parameters: 1009 + linesearch - linesearch context 1010 . steptol - The minimum steplength 1011 . maxstep - The maximum steplength 1012 . rtol - The relative tolerance for iterative line searches 1013 . atol - The absolute tolerance for iterative line searches 1014 . ltol - The change in lambda tolerance for iterative line searches 1015 - max_it - The maximum number of iterations of the line search 1016 1017 Notes: 1018 The user may choose to not set any of the tolerances using PETSC_DEFAULT in 1019 place of an argument. 1020 1021 Level: intermediate 1022 1023 .seealso: SNESLineSearchGetTolerances() 1024 @*/ 1025 PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 1026 { 1027 PetscFunctionBegin; 1028 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1029 PetscValidLogicalCollectiveReal(linesearch,steptol,2); 1030 PetscValidLogicalCollectiveReal(linesearch,maxstep,3); 1031 PetscValidLogicalCollectiveReal(linesearch,rtol,4); 1032 PetscValidLogicalCollectiveReal(linesearch,atol,5); 1033 PetscValidLogicalCollectiveReal(linesearch,ltol,6); 1034 PetscValidLogicalCollectiveInt(linesearch,max_its,7); 1035 1036 if (steptol!= PETSC_DEFAULT) { 1037 if (steptol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %14.12e must be non-negative",(double)steptol); 1038 linesearch->steptol = steptol; 1039 } 1040 1041 if (maxstep!= PETSC_DEFAULT) { 1042 if (maxstep < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %14.12e must be non-negative",(double)maxstep); 1043 linesearch->maxstep = maxstep; 1044 } 1045 1046 if (rtol != PETSC_DEFAULT) { 1047 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); 1048 linesearch->rtol = rtol; 1049 } 1050 1051 if (atol != PETSC_DEFAULT) { 1052 if (atol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %14.12e must be non-negative",(double)atol); 1053 linesearch->atol = atol; 1054 } 1055 1056 if (ltol != PETSC_DEFAULT) { 1057 if (ltol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %14.12e must be non-negative",(double)ltol); 1058 linesearch->ltol = ltol; 1059 } 1060 1061 if (max_its != PETSC_DEFAULT) { 1062 if (max_its < 0) SETERRQ1(PetscObjectComm((PetscObject)linesearch),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its); 1063 linesearch->max_its = max_its; 1064 } 1065 PetscFunctionReturn(0); 1066 } 1067 1068 #undef __FUNCT__ 1069 #define __FUNCT__ "SNESLineSearchGetDamping" 1070 /*@ 1071 SNESLineSearchGetDamping - Gets the line search damping parameter. 1072 1073 Input Parameters: 1074 . linesearch - linesearch context 1075 1076 Output Parameters: 1077 . damping - The damping parameter 1078 1079 Level: advanced 1080 1081 .seealso: SNESLineSearchGetStepTolerance(), SNESQN 1082 @*/ 1083 1084 PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping) 1085 { 1086 PetscFunctionBegin; 1087 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1088 PetscValidPointer(damping, 2); 1089 *damping = linesearch->damping; 1090 PetscFunctionReturn(0); 1091 } 1092 1093 #undef __FUNCT__ 1094 #define __FUNCT__ "SNESLineSearchSetDamping" 1095 /*@ 1096 SNESLineSearchSetDamping - Sets the line search damping paramter. 1097 1098 Input Parameters: 1099 . linesearch - linesearch context 1100 . damping - The damping parameter 1101 1102 Level: intermediate 1103 1104 Notes: 1105 The basic line search merely takes the update step scaled by the damping parameter. 1106 The use of the damping parameter in the l2 and cp line searches is much more subtle; 1107 it is used as a starting point in calculating the secant step. However, the eventual 1108 step may be of greater length than the damping parameter. In the bt line search it is 1109 used as the maximum possible step length, as the bt line search only backtracks. 1110 1111 .seealso: SNESLineSearchGetDamping() 1112 @*/ 1113 PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping) 1114 { 1115 PetscFunctionBegin; 1116 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1117 linesearch->damping = damping; 1118 PetscFunctionReturn(0); 1119 } 1120 1121 #undef __FUNCT__ 1122 #define __FUNCT__ "SNESLineSearchGetOrder" 1123 /*@ 1124 SNESLineSearchGetOrder - Gets the line search approximation order. 1125 1126 Input Parameters: 1127 . linesearch - linesearch context 1128 1129 Output Parameters: 1130 . order - The order 1131 1132 Possible Values for order: 1133 + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 1134 . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 1135 - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 1136 1137 Level: intermediate 1138 1139 .seealso: SNESLineSearchSetOrder() 1140 @*/ 1141 1142 PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch,PetscInt *order) 1143 { 1144 PetscFunctionBegin; 1145 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1146 PetscValidPointer(order, 2); 1147 *order = linesearch->order; 1148 PetscFunctionReturn(0); 1149 } 1150 1151 #undef __FUNCT__ 1152 #define __FUNCT__ "SNESLineSearchSetOrder" 1153 /*@ 1154 SNESLineSearchSetOrder - Sets the line search damping paramter. 1155 1156 Input Parameters: 1157 . linesearch - linesearch context 1158 . order - The damping parameter 1159 1160 Level: intermediate 1161 1162 Possible Values for order: 1163 + 1 or SNES_LINESEARCH_ORDER_LINEAR - linear order 1164 . 2 or SNES_LINESEARCH_ORDER_QUADRATIC - quadratic order 1165 - 3 or SNES_LINESEARCH_ORDER_CUBIC - cubic order 1166 1167 Notes: 1168 Variable orders are supported by the following line searches: 1169 + bt - cubic and quadratic 1170 - cp - linear and quadratic 1171 1172 .seealso: SNESLineSearchGetOrder() 1173 @*/ 1174 PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch,PetscInt order) 1175 { 1176 PetscFunctionBegin; 1177 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1178 linesearch->order = order; 1179 PetscFunctionReturn(0); 1180 } 1181 1182 #undef __FUNCT__ 1183 #define __FUNCT__ "SNESLineSearchGetNorms" 1184 /*@ 1185 SNESLineSearchGetNorms - Gets the norms for for X, Y, and F. 1186 1187 Input Parameters: 1188 . linesearch - linesearch context 1189 1190 Output Parameters: 1191 + xnorm - The norm of the current solution 1192 . fnorm - The norm of the current function 1193 - ynorm - The norm of the current update 1194 1195 Notes: 1196 This function is mainly called from SNES implementations. 1197 1198 Level: developer 1199 1200 .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs() 1201 @*/ 1202 PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm) 1203 { 1204 PetscFunctionBegin; 1205 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1206 if (xnorm) *xnorm = linesearch->xnorm; 1207 if (fnorm) *fnorm = linesearch->fnorm; 1208 if (ynorm) *ynorm = linesearch->ynorm; 1209 PetscFunctionReturn(0); 1210 } 1211 1212 #undef __FUNCT__ 1213 #define __FUNCT__ "SNESLineSearchSetNorms" 1214 /*@ 1215 SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 1216 1217 Input Parameters: 1218 + linesearch - linesearch context 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 Level: advanced 1224 1225 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1226 @*/ 1227 PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 1228 { 1229 PetscFunctionBegin; 1230 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1231 linesearch->xnorm = xnorm; 1232 linesearch->fnorm = fnorm; 1233 linesearch->ynorm = ynorm; 1234 PetscFunctionReturn(0); 1235 } 1236 1237 #undef __FUNCT__ 1238 #define __FUNCT__ "SNESLineSearchComputeNorms" 1239 /*@ 1240 SNESLineSearchComputeNorms - Computes the norms of X, F, and Y. 1241 1242 Input Parameters: 1243 . linesearch - linesearch context 1244 1245 Options Database Keys: 1246 . -snes_linesearch_norms - turn norm computation on or off 1247 1248 Level: intermediate 1249 1250 .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms() 1251 @*/ 1252 PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 1253 { 1254 PetscErrorCode ierr; 1255 SNES snes; 1256 1257 PetscFunctionBegin; 1258 if (linesearch->norms) { 1259 if (linesearch->ops->vinorm) { 1260 ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 1261 ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1262 ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1263 ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr); 1264 } else { 1265 ierr = VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 1266 ierr = VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1267 ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1268 ierr = VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 1269 ierr = VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1270 ierr = VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1271 } 1272 } 1273 PetscFunctionReturn(0); 1274 } 1275 1276 #undef __FUNCT__ 1277 #define __FUNCT__ "SNESLineSearchSetComputeNorms" 1278 /*@ 1279 SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 1280 1281 Input Parameters: 1282 + linesearch - linesearch context 1283 - flg - indicates whether or not to compute norms 1284 1285 Options Database Keys: 1286 . -snes_linesearch_norms - turn norm computation on or off 1287 1288 Notes: 1289 This is most relevant to the SNESLINESEARCHBASIC line search type. 1290 1291 Level: intermediate 1292 1293 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC 1294 @*/ 1295 PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 1296 { 1297 PetscFunctionBegin; 1298 linesearch->norms = flg; 1299 PetscFunctionReturn(0); 1300 } 1301 1302 #undef __FUNCT__ 1303 #define __FUNCT__ "SNESLineSearchGetVecs" 1304 /*@ 1305 SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context 1306 1307 Input Parameters: 1308 . linesearch - linesearch context 1309 1310 Output Parameters: 1311 + X - The old solution 1312 . F - The old function 1313 . Y - The search direction 1314 . W - The new solution 1315 - G - The new function 1316 1317 Level: advanced 1318 1319 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1320 @*/ 1321 PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) 1322 { 1323 PetscFunctionBegin; 1324 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1325 if (X) { 1326 PetscValidPointer(X, 2); 1327 *X = linesearch->vec_sol; 1328 } 1329 if (F) { 1330 PetscValidPointer(F, 3); 1331 *F = linesearch->vec_func; 1332 } 1333 if (Y) { 1334 PetscValidPointer(Y, 4); 1335 *Y = linesearch->vec_update; 1336 } 1337 if (W) { 1338 PetscValidPointer(W, 5); 1339 *W = linesearch->vec_sol_new; 1340 } 1341 if (G) { 1342 PetscValidPointer(G, 6); 1343 *G = linesearch->vec_func_new; 1344 } 1345 PetscFunctionReturn(0); 1346 } 1347 1348 #undef __FUNCT__ 1349 #define __FUNCT__ "SNESLineSearchSetVecs" 1350 /*@ 1351 SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context 1352 1353 Input Parameters: 1354 + linesearch - linesearch context 1355 . X - The old solution 1356 . F - The old function 1357 . Y - The search direction 1358 . W - The new solution 1359 - G - The new function 1360 1361 Level: advanced 1362 1363 .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs() 1364 @*/ 1365 PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) 1366 { 1367 PetscFunctionBegin; 1368 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1369 if (X) { 1370 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1371 linesearch->vec_sol = X; 1372 } 1373 if (F) { 1374 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 1375 linesearch->vec_func = F; 1376 } 1377 if (Y) { 1378 PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 1379 linesearch->vec_update = Y; 1380 } 1381 if (W) { 1382 PetscValidHeaderSpecific(W,VEC_CLASSID,5); 1383 linesearch->vec_sol_new = W; 1384 } 1385 if (G) { 1386 PetscValidHeaderSpecific(G,VEC_CLASSID,6); 1387 linesearch->vec_func_new = G; 1388 } 1389 PetscFunctionReturn(0); 1390 } 1391 1392 #undef __FUNCT__ 1393 #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix" 1394 /*@C 1395 SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1396 SNES options in the database. 1397 1398 Logically Collective on SNESLineSearch 1399 1400 Input Parameters: 1401 + snes - the SNES context 1402 - prefix - the prefix to prepend to all option names 1403 1404 Notes: 1405 A hyphen (-) must NOT be given at the beginning of the prefix name. 1406 The first character of all runtime options is AUTOMATICALLY the hyphen. 1407 1408 Level: advanced 1409 1410 .keywords: SNESLineSearch, append, options, prefix, database 1411 1412 .seealso: SNESGetOptionsPrefix() 1413 @*/ 1414 PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[]) 1415 { 1416 PetscErrorCode ierr; 1417 1418 PetscFunctionBegin; 1419 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1420 ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1421 PetscFunctionReturn(0); 1422 } 1423 1424 #undef __FUNCT__ 1425 #define __FUNCT__ "SNESLineSearchGetOptionsPrefix" 1426 /*@C 1427 SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all 1428 SNESLineSearch options in the database. 1429 1430 Not Collective 1431 1432 Input Parameter: 1433 . linesearch - the SNESLineSearch context 1434 1435 Output Parameter: 1436 . prefix - pointer to the prefix string used 1437 1438 Notes: 1439 On the fortran side, the user should pass in a string 'prefix' of 1440 sufficient length to hold the prefix. 1441 1442 Level: advanced 1443 1444 .keywords: SNESLineSearch, get, options, prefix, database 1445 1446 .seealso: SNESAppendOptionsPrefix() 1447 @*/ 1448 PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[]) 1449 { 1450 PetscErrorCode ierr; 1451 1452 PetscFunctionBegin; 1453 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1454 ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1455 PetscFunctionReturn(0); 1456 } 1457 1458 #undef __FUNCT__ 1459 #define __FUNCT__ "SNESLineSearchSetWorkVecs" 1460 /*@C 1461 SNESLineSearchSetWorkVecs - Gets work vectors for the line search. 1462 1463 Input Parameter: 1464 + linesearch - the SNESLineSearch context 1465 - nwork - the number of work vectors 1466 1467 Level: developer 1468 1469 Developers Note: This is PETSC_EXTERN because it may be used by user written plugin SNES implementations 1470 1471 .keywords: SNESLineSearch, work, vector 1472 1473 .seealso: SNESSetWorkVecs() 1474 @*/ 1475 PetscErrorCode SNESLineSearchSetWorkVecs(SNESLineSearch linesearch, PetscInt nwork) 1476 { 1477 PetscErrorCode ierr; 1478 1479 PetscFunctionBegin; 1480 if (linesearch->vec_sol) { 1481 ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr); 1482 } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 1483 PetscFunctionReturn(0); 1484 } 1485 1486 #undef __FUNCT__ 1487 #define __FUNCT__ "SNESLineSearchGetSuccess" 1488 /*@ 1489 SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application 1490 1491 Input Parameters: 1492 . linesearch - linesearch context 1493 1494 Output Parameters: 1495 . success - The success or failure status 1496 1497 Notes: 1498 This is typically called after SNESLineSearchApply() in order to determine if the line-search failed 1499 (and set the SNES convergence accordingly). 1500 1501 Level: intermediate 1502 1503 .seealso: SNESLineSearchSetSuccess() 1504 @*/ 1505 PetscErrorCode SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success) 1506 { 1507 PetscFunctionBegin; 1508 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1509 PetscValidPointer(success, 2); 1510 if (success) *success = linesearch->success; 1511 PetscFunctionReturn(0); 1512 } 1513 1514 #undef __FUNCT__ 1515 #define __FUNCT__ "SNESLineSearchSetSuccess" 1516 /*@ 1517 SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application 1518 1519 Input Parameters: 1520 + linesearch - linesearch context 1521 - success - The success or failure status 1522 1523 Notes: 1524 This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set 1525 the success or failure of the line search method. 1526 1527 Level: developer 1528 1529 .seealso: SNESLineSearchGetSuccess() 1530 @*/ 1531 PetscErrorCode SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success) 1532 { 1533 PetscFunctionBegin; 1534 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1535 linesearch->success = success; 1536 PetscFunctionReturn(0); 1537 } 1538 1539 #undef __FUNCT__ 1540 #define __FUNCT__ "SNESLineSearchSetVIFunctions" 1541 /*@C 1542 SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 1543 1544 Input Parameters: 1545 + snes - nonlinear context obtained from SNESCreate() 1546 . projectfunc - function for projecting the function to the bounds 1547 - normfunc - function for computing the norm of an active set 1548 1549 Logically Collective on SNES 1550 1551 Calling sequence of projectfunc: 1552 .vb 1553 projectfunc (SNES snes, Vec X) 1554 .ve 1555 1556 Input parameters for projectfunc: 1557 + snes - nonlinear context 1558 - X - current solution 1559 1560 Output parameters for projectfunc: 1561 . X - Projected solution 1562 1563 Calling sequence of normfunc: 1564 .vb 1565 projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm) 1566 .ve 1567 1568 Input parameters for normfunc: 1569 + snes - nonlinear context 1570 . X - current solution 1571 - F - current residual 1572 1573 Output parameters for normfunc: 1574 . fnorm - VI-specific norm of the function 1575 1576 Notes: 1577 The VI solvers require projection of the solution to the feasible set. projectfunc should implement this. 1578 1579 The VI solvers require special evaluation of the function norm such that the norm is only calculated 1580 on the inactive set. This should be implemented by normfunc. 1581 1582 Level: developer 1583 1584 .keywords: SNES, line search, VI, nonlinear, set, line search 1585 1586 .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() 1587 @*/ 1588 extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 1589 { 1590 PetscFunctionBegin; 1591 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1592 if (projectfunc) linesearch->ops->viproject = projectfunc; 1593 if (normfunc) linesearch->ops->vinorm = normfunc; 1594 PetscFunctionReturn(0); 1595 } 1596 1597 #undef __FUNCT__ 1598 #define __FUNCT__ "SNESLineSearchGetVIFunctions" 1599 /*@C 1600 SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 1601 1602 Input Parameters: 1603 . snes - nonlinear context obtained from SNESCreate() 1604 1605 Output Parameters: 1606 + projectfunc - function for projecting the function to the bounds 1607 - normfunc - function for computing the norm of an active set 1608 1609 Logically Collective on SNES 1610 1611 Level: developer 1612 1613 .keywords: SNES, line search, VI, nonlinear, get, line search 1614 1615 .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck() 1616 @*/ 1617 extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 1618 { 1619 PetscFunctionBegin; 1620 if (projectfunc) *projectfunc = linesearch->ops->viproject; 1621 if (normfunc) *normfunc = linesearch->ops->vinorm; 1622 PetscFunctionReturn(0); 1623 } 1624 1625 #undef __FUNCT__ 1626 #define __FUNCT__ "SNESLineSearchRegister" 1627 /*@C 1628 SNESLineSearchRegister - See SNESLineSearchRegisterDynamic() 1629 1630 Level: advanced 1631 @*/ 1632 PetscErrorCode SNESLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNESLineSearch)) 1633 { 1634 char fullname[PETSC_MAX_PATH_LEN]; 1635 PetscErrorCode ierr; 1636 1637 PetscFunctionBegin; 1638 ierr = PetscFunctionListConcat(path,name,fullname);CHKERRQ(ierr); 1639 ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&SNESLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 1640 PetscFunctionReturn(0); 1641 } 1642