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