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