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