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