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_monitor - Print progress of line searches 616 . -snes_linesearch_damping - The linesearch damping parameter 617 . -snes_linesearch_norms - Turn on/off the linesearch norms 618 . -snes_linesearch_keeplambda - Keep the previous search length as the initial guess. 619 . -snes_linesearch_order - The polynomial order of approximation used in the linesearch 620 - -snes_linesearch_max_it - The number of iterations for iterative line searches 621 622 Logically Collective on SNESLineSearch 623 624 Level: intermediate 625 626 627 .seealso: SNESLineSearchCreate() 628 @*/ 629 PetscErrorCode SNESLineSearchSetFromOptions(SNESLineSearch linesearch) { 630 PetscErrorCode ierr; 631 const char *orders[] = {"linear", "quadratic", "cubic"}; 632 PetscInt indx = 0; 633 const char *deft = SNESLINESEARCHBASIC; 634 char type[256]; 635 PetscBool flg, set; 636 PetscFunctionBegin; 637 if (!SNESLineSearchRegisterAllCalled) {ierr = SNESLineSearchRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 638 639 ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr); 640 if (((PetscObject)linesearch)->type_name) { 641 deft = ((PetscObject)linesearch)->type_name; 642 } 643 ierr = PetscOptionsList("-snes_linesearch_type","Line-search method","SNESLineSearchSetType",SNESLineSearchList,deft,type,256,&flg);CHKERRQ(ierr); 644 if (flg) { 645 ierr = SNESLineSearchSetType(linesearch,type);CHKERRQ(ierr); 646 } else if (!((PetscObject)linesearch)->type_name) { 647 ierr = SNESLineSearchSetType(linesearch,deft);CHKERRQ(ierr); 648 } 649 if (linesearch->ops->setfromoptions) { 650 (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr); 651 } 652 653 ierr = PetscOptionsBool("-snes_linesearch_monitor","Print progress of line searches","SNESSNESLineSearchSetMonitor", 654 linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 655 if (set) {ierr = SNESLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);} 656 657 ierr = PetscOptionsReal("-snes_linesearch_damping","Line search damping and initial step guess","SNESLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr); 658 ierr = PetscOptionsReal("-snes_linesearch_rtol","Tolerance for iterative line search","SNESLineSearchSetTolerances",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr); 659 ierr = PetscOptionsBool("-snes_linesearch_norms","Compute final norms in line search","SNESLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr); 660 ierr = PetscOptionsBool("-snes_linesearch_keeplambda","Use previous lambda as damping","SNESLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr); 661 ierr = PetscOptionsInt("-snes_linesearch_max_it","Maximum iterations for iterative line searches","SNESLineSearchSetTolerances",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr); 662 663 ierr = PetscOptionsBool("-snes_linesearch_precheck_picard","Use a correction that sometimes improves convergence of Picard iteration","SNESLineSearchPreCheckPicard",flg,&flg,&set);CHKERRQ(ierr); 664 if (set) { 665 if (flg) { 666 linesearch->precheck_picard_angle = 10.; /* correction only active if angle is less than 10 degrees */ 667 ierr = PetscOptionsReal("-snes_linesearch_precheck_picard_angle","Maximum angle at which to activate the correction", 668 "none",linesearch->precheck_picard_angle,&linesearch->precheck_picard_angle,PETSC_NULL);CHKERRQ(ierr); 669 ierr = SNESLineSearchSetPreCheck(linesearch,SNESLineSearchPreCheckPicard,&linesearch->precheck_picard_angle);CHKERRQ(ierr); 670 } else { 671 ierr = SNESLineSearchSetPreCheck(linesearch,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 672 } 673 } 674 ierr = PetscOptionsEList("-snes_linesearch_order", "Order of approximation", "SNESLineSearchSetOrder", orders,3,orders[(PetscInt)linesearch->order],&indx,&flg);CHKERRQ(ierr); 675 if (flg) { 676 switch (indx) { 677 case 0: linesearch->order = SNES_LINESEARCH_LINEAR; 678 break; 679 case 1: linesearch->order = SNES_LINESEARCH_QUADRATIC; 680 break; 681 case 2: linesearch->order = SNES_LINESEARCH_CUBIC; 682 break; 683 } 684 } 685 686 687 ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr); 688 ierr = PetscOptionsEnd();CHKERRQ(ierr); 689 PetscFunctionReturn(0); 690 } 691 692 #undef __FUNCT__ 693 #define __FUNCT__ "SNESLineSearchView" 694 /*@ 695 SNESLineSearchView - Prints useful information about the line search not 696 related to an individual call. 697 698 Input Parameters: 699 . linesearch - linesearch context 700 701 Logically Collective on SNESLineSearch 702 703 Level: intermediate 704 705 .seealso: SNESLineSearchCreate() 706 @*/ 707 PetscErrorCode SNESLineSearchView(SNESLineSearch linesearch, PetscViewer viewer) { 708 PetscErrorCode ierr; 709 PetscBool iascii; 710 PetscFunctionBegin; 711 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 712 if (!viewer) { 713 ierr = PetscViewerASCIIGetStdout(((PetscObject)linesearch)->comm,&viewer);CHKERRQ(ierr); 714 } 715 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 716 PetscCheckSameComm(linesearch,1,viewer,2); 717 718 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 719 if (iascii) { 720 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)linesearch,viewer,"SNESLineSearch Object");CHKERRQ(ierr); 721 if (linesearch->ops->view) { 722 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 723 ierr = (*linesearch->ops->view)(linesearch,viewer);CHKERRQ(ierr); 724 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 725 } 726 ierr = PetscViewerASCIIPrintf(viewer," maxstep=%G, minlambda=%G\n", linesearch->maxstep,linesearch->steptol);CHKERRQ(ierr); 727 ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%G, absolute=%G, lambda=%G\n", linesearch->rtol,linesearch->atol,linesearch->ltol);CHKERRQ(ierr); 728 ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D\n", linesearch->max_its);CHKERRQ(ierr); 729 if (linesearch->ops->precheckstep) { 730 if (linesearch->ops->precheckstep == SNESLineSearchPreCheckPicard) { 731 ierr = PetscViewerASCIIPrintf(viewer," using precheck step to speed up Picard convergence\n", linesearch->max_its);CHKERRQ(ierr); 732 } else { 733 ierr = PetscViewerASCIIPrintf(viewer," using user-defined precheck step\n", linesearch->max_its);CHKERRQ(ierr); 734 } 735 } 736 if (linesearch->ops->postcheckstep) { 737 ierr = PetscViewerASCIIPrintf(viewer," using user-defined postcheck step\n", linesearch->max_its);CHKERRQ(ierr); 738 } 739 } 740 PetscFunctionReturn(0); 741 } 742 743 #undef __FUNCT__ 744 #define __FUNCT__ "SNESLineSearchSetType" 745 /*@C 746 SNESLineSearchSetType - Sets the linesearch type 747 748 Input Parameters: 749 + linesearch - linesearch context 750 - type - The type of line search to be used 751 752 Available Types: 753 + basic - Simple damping line search. 754 . bt - Backtracking line search over the L2 norm of the function 755 . l2 - Secant line search over the L2 norm of the function 756 . cp - Critical point secant line search assuming F(x) = grad G(x) for some unknown G(x) 757 - shell - User provided SNESLineSearch implementation 758 759 Logically Collective on SNESLineSearch 760 761 Level: intermediate 762 763 764 .seealso: SNESLineSearchCreate() 765 @*/ 766 PetscErrorCode SNESLineSearchSetType(SNESLineSearch linesearch, const SNESLineSearchType type) 767 { 768 769 PetscErrorCode ierr,(*r)(SNESLineSearch); 770 PetscBool match; 771 772 PetscFunctionBegin; 773 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 774 PetscValidCharPointer(type,2); 775 776 ierr = PetscTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr); 777 if (match) PetscFunctionReturn(0); 778 779 ierr = PetscFListFind(SNESLineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 780 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type); 781 /* Destroy the previous private linesearch context */ 782 if (linesearch->ops->destroy) { 783 ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr); 784 linesearch->ops->destroy = PETSC_NULL; 785 } 786 /* Reinitialize function pointers in SNESLineSearchOps structure */ 787 linesearch->ops->apply = 0; 788 linesearch->ops->view = 0; 789 linesearch->ops->setfromoptions = 0; 790 linesearch->ops->destroy = 0; 791 792 ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr); 793 ierr = (*r)(linesearch);CHKERRQ(ierr); 794 #if defined(PETSC_HAVE_AMS) 795 if (PetscAMSPublishAll) { 796 ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr); 797 } 798 #endif 799 PetscFunctionReturn(0); 800 } 801 802 #undef __FUNCT__ 803 #define __FUNCT__ "SNESLineSearchSetSNES" 804 /*@ 805 SNESLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation. 806 807 Input Parameters: 808 + linesearch - linesearch context 809 - snes - The snes instance 810 811 Level: developer 812 813 Notes: 814 This happens automatically when the line search is gotten/created with 815 SNESGetSNESLineSearch(). This routine is therefore mainly called within SNES 816 implementations. 817 818 819 .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 820 @*/ 821 PetscErrorCode SNESLineSearchSetSNES(SNESLineSearch linesearch, SNES snes){ 822 PetscFunctionBegin; 823 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 824 PetscValidHeaderSpecific(snes,SNES_CLASSID,2); 825 linesearch->snes = snes; 826 PetscFunctionReturn(0); 827 } 828 829 #undef __FUNCT__ 830 #define __FUNCT__ "SNESLineSearchGetSNES" 831 /*@ 832 SNESLineSearchGetSNES - Gets the SNES for the linesearch for function evaluation. 833 834 Input Parameters: 835 . linesearch - linesearch context 836 837 Output Parameters: 838 . snes - The snes instance 839 840 Level: advanced 841 842 .seealso: SNESLineSearchGetSNES(), SNESLineSearchSetVecs(), SNES 843 @*/ 844 PetscErrorCode SNESLineSearchGetSNES(SNESLineSearch linesearch, SNES *snes){ 845 PetscFunctionBegin; 846 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 847 PetscValidPointer(snes, 2); 848 *snes = linesearch->snes; 849 PetscFunctionReturn(0); 850 } 851 852 #undef __FUNCT__ 853 #define __FUNCT__ "SNESLineSearchGetLambda" 854 /*@ 855 SNESLineSearchGetLambda - Gets the last linesearch steplength discovered. 856 857 Input Parameters: 858 . linesearch - linesearch context 859 860 Output Parameters: 861 . lambda - The last steplength computed during SNESLineSearchApply() 862 863 Level: advanced 864 865 Notes: 866 This is useful in methods where the solver is ill-scaled and 867 requires some adaptive notion of the difference in scale between the 868 solution and the function. For instance, SNESQN may be scaled by the 869 line search lambda using the argument -snes_qn_scaling ls. 870 871 872 .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetDamping(), SNESLineSearchApply() 873 @*/ 874 PetscErrorCode SNESLineSearchGetLambda(SNESLineSearch linesearch,PetscReal *lambda) 875 { 876 PetscFunctionBegin; 877 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 878 PetscValidPointer(lambda, 2); 879 *lambda = linesearch->lambda; 880 PetscFunctionReturn(0); 881 } 882 883 #undef __FUNCT__ 884 #define __FUNCT__ "SNESLineSearchSetLambda" 885 /*@ 886 SNESLineSearchSetLambda - Sets the linesearch steplength. 887 888 Input Parameters: 889 + linesearch - linesearch context 890 - lambda - The last steplength. 891 892 Notes: 893 This routine is typically used within implementations of SNESLineSearchApply 894 to set the final steplength. This routine (and SNESLineSearchGetLambda()) were 895 added in order to facilitate Quasi-Newton methods that use the previous steplength 896 as an inner scaling parameter. 897 898 Level: advanced 899 900 .seealso: SNESLineSearchGetLambda() 901 @*/ 902 PetscErrorCode SNESLineSearchSetLambda(SNESLineSearch linesearch, PetscReal lambda) 903 { 904 PetscFunctionBegin; 905 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 906 linesearch->lambda = lambda; 907 PetscFunctionReturn(0); 908 } 909 910 #undef __FUNCT__ 911 #define __FUNCT__ "SNESLineSearchGetTolerances" 912 /*@ 913 SNESLineSearchGetTolerances - Gets the tolerances for the method. These include 914 tolerances for the relative and absolute change in the function norm, the change 915 in lambda for iterative line searches, the minimum steplength, the maximum steplength, 916 and the maximum number of iterations the line search procedure may take. 917 918 Input Parameters: 919 . linesearch - linesearch context 920 921 Output Parameters: 922 + steptol - The minimum steplength 923 . maxstep - The maximum steplength 924 . rtol - The relative tolerance for iterative line searches 925 . atol - The absolute tolerance for iterative line searches 926 . ltol - The change in lambda tolerance for iterative line searches 927 - max_it - The maximum number of iterations of the line search 928 929 Level: intermediate 930 931 Notes: 932 Different line searches may implement these parameters slightly differently as 933 the method requires. 934 935 .seealso: SNESLineSearchSetTolerances() 936 @*/ 937 PetscErrorCode SNESLineSearchGetTolerances(SNESLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 938 { 939 PetscFunctionBegin; 940 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 941 if (steptol) { 942 PetscValidPointer(steptol, 2); 943 *steptol = linesearch->steptol; 944 } 945 if (maxstep) { 946 PetscValidPointer(maxstep, 3); 947 *maxstep = linesearch->maxstep; 948 } 949 if (rtol) { 950 PetscValidPointer(rtol, 4); 951 *rtol = linesearch->rtol; 952 } 953 if (atol) { 954 PetscValidPointer(atol, 5); 955 *atol = linesearch->atol; 956 } 957 if (ltol) { 958 PetscValidPointer(ltol, 6); 959 *ltol = linesearch->ltol; 960 } 961 if (max_its) { 962 PetscValidPointer(max_its, 7); 963 *max_its = linesearch->max_its; 964 } 965 PetscFunctionReturn(0); 966 } 967 968 #undef __FUNCT__ 969 #define __FUNCT__ "SNESLineSearchSetTolerances" 970 /*@ 971 SNESLineSearchSetTolerances - Gets the tolerances for the method. These include 972 tolerances for the relative and absolute change in the function norm, the change 973 in lambda for iterative line searches, the minimum steplength, the maximum steplength, 974 and the maximum number of iterations the line search procedure may take. 975 976 Input Parameters: 977 + linesearch - linesearch context 978 . steptol - The minimum steplength 979 . maxstep - The maximum steplength 980 . rtol - The relative tolerance for iterative line searches 981 . atol - The absolute tolerance for iterative line searches 982 . ltol - The change in lambda tolerance for iterative line searches 983 - max_it - The maximum number of iterations of the line search 984 985 Notes: 986 The user may choose to not set any of the tolerances in the method using PETSC_DEFAULT in 987 place of an argument. 988 989 Level: intermediate 990 991 .seealso: SNESLineSearchGetTolerances() 992 @*/ 993 PetscErrorCode SNESLineSearchSetTolerances(SNESLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 994 { 995 PetscFunctionBegin; 996 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 997 PetscValidLogicalCollectiveReal(linesearch,steptol,2); 998 PetscValidLogicalCollectiveReal(linesearch,maxstep,3); 999 PetscValidLogicalCollectiveReal(linesearch,rtol,4); 1000 PetscValidLogicalCollectiveReal(linesearch,atol,5); 1001 PetscValidLogicalCollectiveReal(linesearch,ltol,6); 1002 PetscValidLogicalCollectiveInt(linesearch,max_its,7); 1003 1004 if ( steptol!= PETSC_DEFAULT) { 1005 if ( steptol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Minimum step length %G must be non-negative",steptol); 1006 linesearch->steptol = steptol; 1007 } 1008 1009 if ( maxstep!= PETSC_DEFAULT) { 1010 if ( maxstep < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum step length %G must be non-negative",maxstep); 1011 linesearch->maxstep = maxstep; 1012 } 1013 1014 if (rtol != PETSC_DEFAULT) { 1015 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); 1016 linesearch->rtol = rtol; 1017 } 1018 1019 if (atol != PETSC_DEFAULT) { 1020 if (atol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %G must be non-negative",atol); 1021 linesearch->atol = atol; 1022 } 1023 1024 if (ltol != PETSC_DEFAULT) { 1025 if (ltol < 0.0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Labmda tolerance %G must be non-negative",ltol); 1026 linesearch->ltol = ltol; 1027 } 1028 1029 if (max_its != PETSC_DEFAULT) { 1030 if (max_its < 0) SETERRQ1(((PetscObject)linesearch)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",max_its); 1031 linesearch->max_its = max_its; 1032 } 1033 1034 PetscFunctionReturn(0); 1035 } 1036 1037 1038 #undef __FUNCT__ 1039 #define __FUNCT__ "SNESLineSearchGetDamping" 1040 /*@ 1041 SNESLineSearchGetDamping - Gets the line search damping parameter. 1042 1043 Input Parameters: 1044 . linesearch - linesearch context 1045 1046 Output Parameters: 1047 . damping - The damping parameter 1048 1049 Level: advanced 1050 1051 .seealso: SNESLineSearchGetStepTolerance(), SNESQN 1052 @*/ 1053 1054 PetscErrorCode SNESLineSearchGetDamping(SNESLineSearch linesearch,PetscReal *damping) 1055 { 1056 PetscFunctionBegin; 1057 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1058 PetscValidPointer(damping, 2); 1059 *damping = linesearch->damping; 1060 PetscFunctionReturn(0); 1061 } 1062 1063 #undef __FUNCT__ 1064 #define __FUNCT__ "SNESLineSearchSetDamping" 1065 /*@ 1066 SNESLineSearchSetDamping - Sets the line search damping paramter. 1067 1068 Input Parameters: 1069 . linesearch - linesearch context 1070 . damping - The damping parameter 1071 1072 Level: intermediate 1073 1074 Notes: 1075 The basic line search merely takes the update step scaled by the damping parameter. 1076 The use of the damping parameter in the l2 and cp line searches is much more subtle; 1077 it is used as a starting point in calculating the secant step. However, the eventual 1078 step may be of greater length than the damping parameter. In the bt line search it is 1079 used as the maximum possible step length, as the bt line search only backtracks. 1080 1081 .seealso: SNESLineSearchGetDamping() 1082 @*/ 1083 PetscErrorCode SNESLineSearchSetDamping(SNESLineSearch linesearch,PetscReal damping) 1084 { 1085 PetscFunctionBegin; 1086 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1087 linesearch->damping = damping; 1088 PetscFunctionReturn(0); 1089 } 1090 1091 #undef __FUNCT__ 1092 #define __FUNCT__ "SNESLineSearchGetOrder" 1093 /*@ 1094 SNESLineSearchGetOrder - Gets the line search approximation order. 1095 1096 Input Parameters: 1097 . linesearch - linesearch context 1098 1099 Output Parameters: 1100 . order - The order 1101 1102 Possible Values for order: 1103 + SNES_LINESEARCH_LINEAR - linear method 1104 . SNES_LINESEARCH_QUADRATIC - quadratic method 1105 - SNES_LINESEARCH_CUBIC - cubic method 1106 1107 Level: intermediate 1108 1109 .seealso: SNESLineSearchSetOrder() 1110 @*/ 1111 1112 PetscErrorCode SNESLineSearchGetOrder(SNESLineSearch linesearch,SNESLineSearchOrder *order) 1113 { 1114 PetscFunctionBegin; 1115 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1116 PetscValidPointer(order, 2); 1117 *order = linesearch->order; 1118 PetscFunctionReturn(0); 1119 } 1120 1121 #undef __FUNCT__ 1122 #define __FUNCT__ "SNESLineSearchSetOrder" 1123 /*@ 1124 SNESLineSearchSetOrder - Sets the line search damping paramter. 1125 1126 Input Parameters: 1127 . linesearch - linesearch context 1128 . order - The damping parameter 1129 1130 Level: intermediate 1131 1132 Possible Values for order: 1133 + SNES_LINESEARCH_LINEAR - linear method 1134 . SNES_LINESEARCH_QUADRATIC - quadratic method 1135 - SNES_LINESEARCH_CUBIC - cubic method 1136 1137 Notes: 1138 Variable orders are supported by the following line searches: 1139 + bt - cubic and quadratic 1140 - cp - linear and quadratic 1141 1142 .seealso: SNESLineSearchGetOrder() 1143 @*/ 1144 PetscErrorCode SNESLineSearchSetOrder(SNESLineSearch linesearch,SNESLineSearchOrder order) 1145 { 1146 PetscFunctionBegin; 1147 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1148 linesearch->order = order; 1149 PetscFunctionReturn(0); 1150 } 1151 1152 #undef __FUNCT__ 1153 #define __FUNCT__ "SNESLineSearchGetNorms" 1154 /*@ 1155 SNESLineSearchGetNorms - Gets the norms for for X, Y, and F. 1156 1157 Input Parameters: 1158 . linesearch - linesearch context 1159 1160 Output Parameters: 1161 + xnorm - The norm of the current solution 1162 . fnorm - The norm of the current function 1163 - ynorm - The norm of the current update 1164 1165 Notes: 1166 This function is mainly called from SNES implementations. 1167 1168 Level: developer 1169 1170 .seealso: SNESLineSearchSetNorms() SNESLineSearchGetVecs() 1171 @*/ 1172 PetscErrorCode SNESLineSearchGetNorms(SNESLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm) 1173 { 1174 PetscFunctionBegin; 1175 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1176 if (xnorm) { 1177 *xnorm = linesearch->xnorm; 1178 } 1179 if (fnorm) { 1180 *fnorm = linesearch->fnorm; 1181 } 1182 if (ynorm) { 1183 *ynorm = linesearch->ynorm; 1184 } 1185 PetscFunctionReturn(0); 1186 } 1187 1188 #undef __FUNCT__ 1189 #define __FUNCT__ "SNESLineSearchSetNorms" 1190 /*@ 1191 SNESLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 1192 1193 Input Parameters: 1194 + linesearch - linesearch context 1195 . xnorm - The norm of the current solution 1196 . fnorm - The norm of the current function 1197 - ynorm - The norm of the current update 1198 1199 Level: advanced 1200 1201 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1202 @*/ 1203 PetscErrorCode SNESLineSearchSetNorms(SNESLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 1204 { 1205 PetscFunctionBegin; 1206 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1207 linesearch->xnorm = xnorm; 1208 linesearch->fnorm = fnorm; 1209 linesearch->ynorm = ynorm; 1210 PetscFunctionReturn(0); 1211 } 1212 1213 #undef __FUNCT__ 1214 #define __FUNCT__ "SNESLineSearchComputeNorms" 1215 /*@ 1216 SNESLineSearchComputeNorms - Computes the norms of X, F, and Y. 1217 1218 Input Parameters: 1219 . linesearch - linesearch context 1220 1221 Options Database Keys: 1222 . -snes_linesearch_norms - turn norm computation on or off 1223 1224 Level: intermediate 1225 1226 .seealso: SNESLineSearchGetNorms, SNESLineSearchSetNorms(), SNESLineSearchSetComputeNorms() 1227 @*/ 1228 PetscErrorCode SNESLineSearchComputeNorms(SNESLineSearch linesearch) 1229 { 1230 PetscErrorCode ierr; 1231 SNES snes; 1232 PetscFunctionBegin; 1233 if (linesearch->norms) { 1234 if (linesearch->ops->vinorm) { 1235 ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 1236 ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1237 ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1238 ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr); 1239 } else { 1240 ierr = VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 1241 ierr = VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1242 ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1243 ierr = VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 1244 ierr = VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1245 ierr = VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1246 } 1247 } 1248 PetscFunctionReturn(0); 1249 } 1250 1251 1252 #undef __FUNCT__ 1253 #define __FUNCT__ "SNESLineSearchSetComputeNorms" 1254 /*@ 1255 SNESLineSearchSetComputeNorms - Turns on or off the computation of final norms in the line search. 1256 1257 Input Parameters: 1258 + linesearch - linesearch context 1259 - flg - indicates whether or not to compute norms 1260 1261 Options Database Keys: 1262 . -snes_linesearch_norms - turn norm computation on or off 1263 1264 Notes: 1265 This is most relevant to the SNESLINESEARCHBASIC line search type. 1266 1267 Level: intermediate 1268 1269 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetNorms(), SNESLineSearchComputeNorms(), SNESLINESEARCHBASIC 1270 @*/ 1271 PetscErrorCode SNESLineSearchSetComputeNorms(SNESLineSearch linesearch, PetscBool flg) 1272 { 1273 PetscFunctionBegin; 1274 linesearch->norms = flg; 1275 PetscFunctionReturn(0); 1276 } 1277 1278 #undef __FUNCT__ 1279 #define __FUNCT__ "SNESLineSearchGetVecs" 1280 /*@ 1281 SNESLineSearchGetVecs - Gets the vectors from the SNESLineSearch context 1282 1283 Input Parameters: 1284 . linesearch - linesearch context 1285 1286 Output Parameters: 1287 + X - The old solution 1288 . F - The old function 1289 . Y - The search direction 1290 . W - The new solution 1291 - G - The new function 1292 1293 Level: advanced 1294 1295 .seealso: SNESLineSearchGetNorms(), SNESLineSearchSetVecs() 1296 @*/ 1297 PetscErrorCode SNESLineSearchGetVecs(SNESLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) { 1298 PetscFunctionBegin; 1299 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1300 if (X) { 1301 PetscValidPointer(X, 2); 1302 *X = linesearch->vec_sol; 1303 } 1304 if (F) { 1305 PetscValidPointer(F, 3); 1306 *F = linesearch->vec_func; 1307 } 1308 if (Y) { 1309 PetscValidPointer(Y, 4); 1310 *Y = linesearch->vec_update; 1311 } 1312 if (W) { 1313 PetscValidPointer(W, 5); 1314 *W = linesearch->vec_sol_new; 1315 } 1316 if (G) { 1317 PetscValidPointer(G, 6); 1318 *G = linesearch->vec_func_new; 1319 } 1320 1321 PetscFunctionReturn(0); 1322 } 1323 1324 #undef __FUNCT__ 1325 #define __FUNCT__ "SNESLineSearchSetVecs" 1326 /*@ 1327 SNESLineSearchSetVecs - Sets the vectors on the SNESLineSearch context 1328 1329 Input Parameters: 1330 + linesearch - linesearch context 1331 . X - The old solution 1332 . F - The old function 1333 . Y - The search direction 1334 . W - The new solution 1335 - G - The new function 1336 1337 Level: advanced 1338 1339 .seealso: SNESLineSearchSetNorms(), SNESLineSearchGetVecs() 1340 @*/ 1341 PetscErrorCode SNESLineSearchSetVecs(SNESLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) { 1342 PetscFunctionBegin; 1343 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1344 if (X) { 1345 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1346 linesearch->vec_sol = X; 1347 } 1348 if (F) { 1349 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 1350 linesearch->vec_func = F; 1351 } 1352 if (Y) { 1353 PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 1354 linesearch->vec_update = Y; 1355 } 1356 if (W) { 1357 PetscValidHeaderSpecific(W,VEC_CLASSID,5); 1358 linesearch->vec_sol_new = W; 1359 } 1360 if (G) { 1361 PetscValidHeaderSpecific(G,VEC_CLASSID,6); 1362 linesearch->vec_func_new = G; 1363 } 1364 1365 PetscFunctionReturn(0); 1366 } 1367 1368 #undef __FUNCT__ 1369 #define __FUNCT__ "SNESLineSearchAppendOptionsPrefix" 1370 /*@C 1371 SNESLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1372 SNES options in the database. 1373 1374 Logically Collective on SNESLineSearch 1375 1376 Input Parameters: 1377 + snes - the SNES context 1378 - prefix - the prefix to prepend to all option names 1379 1380 Notes: 1381 A hyphen (-) must NOT be given at the beginning of the prefix name. 1382 The first character of all runtime options is AUTOMATICALLY the hyphen. 1383 1384 Level: advanced 1385 1386 .keywords: SNESLineSearch, append, options, prefix, database 1387 1388 .seealso: SNESGetOptionsPrefix() 1389 @*/ 1390 PetscErrorCode SNESLineSearchAppendOptionsPrefix(SNESLineSearch linesearch,const char prefix[]) 1391 { 1392 PetscErrorCode ierr; 1393 1394 PetscFunctionBegin; 1395 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1396 ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1397 PetscFunctionReturn(0); 1398 } 1399 1400 #undef __FUNCT__ 1401 #define __FUNCT__ "SNESLineSearchGetOptionsPrefix" 1402 /*@C 1403 SNESLineSearchGetOptionsPrefix - Sets the prefix used for searching for all 1404 SNESLineSearch options in the database. 1405 1406 Not Collective 1407 1408 Input Parameter: 1409 . linesearch - the SNESLineSearch context 1410 1411 Output Parameter: 1412 . prefix - pointer to the prefix string used 1413 1414 Notes: 1415 On the fortran side, the user should pass in a string 'prefix' of 1416 sufficient length to hold the prefix. 1417 1418 Level: advanced 1419 1420 .keywords: SNESLineSearch, get, options, prefix, database 1421 1422 .seealso: SNESAppendOptionsPrefix() 1423 @*/ 1424 PetscErrorCode SNESLineSearchGetOptionsPrefix(SNESLineSearch linesearch,const char *prefix[]) 1425 { 1426 PetscErrorCode ierr; 1427 1428 PetscFunctionBegin; 1429 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1430 ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1431 PetscFunctionReturn(0); 1432 } 1433 1434 #undef __FUNCT__ 1435 #define __FUNCT__ "SNESLineSearchGetWork" 1436 /*@ 1437 SNESLineSearchGetWork - Gets work vectors for the line search. 1438 1439 Input Parameter: 1440 + linesearch - the SNESLineSearch context 1441 - nwork - the number of work vectors 1442 1443 Level: developer 1444 1445 Notes: 1446 This is typically called at the beginning of a SNESLineSearch or SNESLineSearchShell implementation. 1447 1448 .keywords: SNESLineSearch, work, vector 1449 1450 .seealso: SNESDefaultGetWork() 1451 @*/ 1452 PetscErrorCode SNESLineSearchGetWork(SNESLineSearch linesearch, PetscInt nwork) 1453 { 1454 PetscErrorCode ierr; 1455 PetscFunctionBegin; 1456 if (linesearch->vec_sol) { 1457 ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr); 1458 } else { 1459 SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 1460 } 1461 PetscFunctionReturn(0); 1462 } 1463 1464 1465 #undef __FUNCT__ 1466 #define __FUNCT__ "SNESLineSearchGetSuccess" 1467 /*@ 1468 SNESLineSearchGetSuccess - Gets the success/failure status of the last line search application 1469 1470 Input Parameters: 1471 . linesearch - linesearch context 1472 1473 Output Parameters: 1474 . success - The success or failure status 1475 1476 Notes: 1477 This is typically called after SNESLineSearchApply in order to determine if the line-search failed 1478 (and set the SNES convergence accordingly). 1479 1480 Level: intermediate 1481 1482 .seealso: SNESLineSearchSetSuccess() 1483 @*/ 1484 PetscErrorCode SNESLineSearchGetSuccess(SNESLineSearch linesearch, PetscBool *success) 1485 { 1486 PetscFunctionBegin; 1487 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1488 PetscValidPointer(success, 2); 1489 if (success) { 1490 *success = linesearch->success; 1491 } 1492 PetscFunctionReturn(0); 1493 } 1494 1495 #undef __FUNCT__ 1496 #define __FUNCT__ "SNESLineSearchSetSuccess" 1497 /*@ 1498 SNESLineSearchSetSuccess - Sets the success/failure status of the last line search application 1499 1500 Input Parameters: 1501 + linesearch - linesearch context 1502 - success - The success or failure status 1503 1504 Notes: 1505 This is typically called in a SNESLineSearchApply() or SNESLineSearchShell implementation to set 1506 the success or failure of the line search method. 1507 1508 Level: developer 1509 1510 .seealso: SNESLineSearchGetSuccess() 1511 @*/ 1512 PetscErrorCode SNESLineSearchSetSuccess(SNESLineSearch linesearch, PetscBool success) 1513 { 1514 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1515 PetscFunctionBegin; 1516 linesearch->success = success; 1517 PetscFunctionReturn(0); 1518 } 1519 1520 #undef __FUNCT__ 1521 #define __FUNCT__ "SNESLineSearchSetVIFunctions" 1522 /*@C 1523 SNESLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 1524 1525 Input Parameters: 1526 + snes - nonlinear context obtained from SNESCreate() 1527 . projectfunc - function for projecting the function to the bounds 1528 - normfunc - function for computing the norm of an active set 1529 1530 Logically Collective on SNES 1531 1532 Calling sequence of projectfunc: 1533 .vb 1534 projectfunc (SNES snes, Vec X) 1535 .ve 1536 1537 Input parameters for projectfunc: 1538 + snes - nonlinear context 1539 - X - current solution 1540 1541 Output parameters for projectfunc: 1542 . X - Projected solution 1543 1544 Calling sequence of normfunc: 1545 .vb 1546 projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm) 1547 .ve 1548 1549 Input parameters for normfunc: 1550 + snes - nonlinear context 1551 . X - current solution 1552 - F - current residual 1553 1554 Output parameters for normfunc: 1555 . fnorm - VI-specific norm of the function 1556 1557 Notes: 1558 The VI solvers require projection of the solution to the feasible set. projectfunc should implement this. 1559 1560 The VI solvers require special evaluation of the function norm such that the norm is only calculated 1561 on the inactive set. This should be implemented by normfunc. 1562 1563 Level: developer 1564 1565 .keywords: SNES, line search, VI, nonlinear, set, line search 1566 1567 .seealso: SNESLineSearchGetVIFunctions(), SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() 1568 @*/ 1569 extern PetscErrorCode SNESLineSearchSetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc projectfunc, SNESLineSearchVINormFunc normfunc) 1570 { 1571 PetscFunctionBegin; 1572 PetscValidHeaderSpecific(linesearch,SNESLINESEARCH_CLASSID,1); 1573 if (projectfunc) linesearch->ops->viproject = projectfunc; 1574 if (normfunc) linesearch->ops->vinorm = normfunc; 1575 PetscFunctionReturn(0); 1576 } 1577 1578 #undef __FUNCT__ 1579 #define __FUNCT__ "SNESLineSearchGetVIFunctions" 1580 /*@C 1581 SNESLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 1582 1583 Input Parameters: 1584 . snes - nonlinear context obtained from SNESCreate() 1585 1586 Output Parameters: 1587 + projectfunc - function for projecting the function to the bounds 1588 - normfunc - function for computing the norm of an active set 1589 1590 Logically Collective on SNES 1591 1592 Level: developer 1593 1594 .keywords: SNES, line search, VI, nonlinear, get, line search 1595 1596 .seealso: SNESLineSearchSetVIFunctions(), SNESLineSearchGetPostCheck(), SNESLineSearchGetPreCheck() 1597 @*/ 1598 extern PetscErrorCode SNESLineSearchGetVIFunctions(SNESLineSearch linesearch, SNESLineSearchVIProjectFunc *projectfunc, SNESLineSearchVINormFunc *normfunc) 1599 { 1600 PetscFunctionBegin; 1601 if (projectfunc) *projectfunc = linesearch->ops->viproject; 1602 if (normfunc) *normfunc = linesearch->ops->vinorm; 1603 PetscFunctionReturn(0); 1604 } 1605 1606 #undef __FUNCT__ 1607 #define __FUNCT__ "SNESLineSearchRegister" 1608 /*@C 1609 SNESLineSearchRegister - See SNESLineSearchRegisterDynamic() 1610 1611 Level: advanced 1612 @*/ 1613 PetscErrorCode SNESLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNESLineSearch)) 1614 { 1615 char fullname[PETSC_MAX_PATH_LEN]; 1616 PetscErrorCode ierr; 1617 1618 PetscFunctionBegin; 1619 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 1620 ierr = PetscFListAdd(&SNESLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 1621 PetscFunctionReturn(0); 1622 } 1623