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