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