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