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