1 #include <private/linesearchimpl.h> /*I "petscsnes.h" I*/ 2 3 PetscBool PetscLineSearchRegisterAllCalled = PETSC_FALSE; 4 PetscFList PetscLineSearchList = PETSC_NULL; 5 6 PetscClassId PETSCLINESEARCH_CLASSID; 7 PetscLogEvent PetscLineSearch_Apply; 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "PetscLineSearchCreate" 11 /*@ 12 PetscLineSearchCreate - 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 PetscLineSearchCreate(MPI_Comm comm, PetscLineSearch *outlinesearch) { 30 PetscErrorCode ierr; 31 PetscLineSearch linesearch; 32 PetscFunctionBegin; 33 PetscValidPointer(outlinesearch,2); 34 *outlinesearch = PETSC_NULL; 35 ierr = PetscHeaderCreate(linesearch,_p_LineSearch,struct _LineSearchOps,PETSCLINESEARCH_CLASSID, 0, 36 "PetscLineSearch","Line-search method","PetscLineSearch",comm,PetscLineSearchDestroy,PetscLineSearchView);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__ "PetscLineSearchSetUp" 70 /*@ 71 PetscLineSearchSetUp - 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: PetscLineSearch, SetUp 81 82 .seealso: PetscLineSearchReset() 83 @*/ 84 85 PetscErrorCode PetscLineSearchSetUp(PetscLineSearch linesearch) { 86 PetscErrorCode ierr; 87 PetscFunctionBegin; 88 if (!((PetscObject)linesearch)->type_name) { 89 ierr = PetscLineSearchSetType(linesearch,PETSCLINESEARCHBASIC);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__ "PetscLineSearchReset" 109 110 /*@ 111 PetscLineSearchReset - Tears down the structures required for application 112 113 Collective on PetscLineSearch 114 115 Input Parameters: 116 . linesearch - The LineSearch instance. 117 118 Level: Intermediate 119 120 .keywords: PetscLineSearch, Create 121 122 .seealso: PetscLineSearchSetUp() 123 @*/ 124 125 PetscErrorCode PetscLineSearchReset(PetscLineSearch 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__ "PetscLineSearchSetPreCheck" 143 /*@C 144 PetscLineSearchSetPreCheck - Sets a pre-check function for the line search routine. 145 146 Logically Collective on PetscLineSearch 147 148 Input Parameters: 149 + linesearch - the PetscLineSearch 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 (PetscLineSearch 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: PetscLineSearchSetPostCheck() 166 @*/ 167 PetscErrorCode PetscLineSearchSetPreCheck(PetscLineSearch linesearch, PetscLineSearchPreCheckFunc func,void *ctx) 168 { 169 PetscFunctionBegin; 170 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_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__ "PetscLineSearchGetPreCheck" 179 /*@C 180 PetscLineSearchSetPreCheck - Sets a pre-check function for the line search routine. 181 182 Input Parameters: 183 . linesearch - the PetscLineSearch 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: PetscLineSearchGetPostCheck(), PetscLineSearchSetPreCheck() 195 @*/ 196 PetscErrorCode PetscLineSearchGetPreCheck(PetscLineSearch linesearch, PetscLineSearchPreCheckFunc *func,void **ctx) 197 { 198 PetscFunctionBegin; 199 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_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__ "PetscLineSearchSetPostCheck" 208 /*@C 209 PetscLineSearchSetPostCheck - Sets a post-check function for the line search routine. 210 211 Logically Collective on PetscLineSearch 212 213 Input Parameters: 214 + linesearch - the PetscLineSearch 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 (PetscLineSearch 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: PetscLineSearchSetPreCheck() 233 @*/ 234 PetscErrorCode PetscLineSearchSetPostCheck(PetscLineSearch linesearch, PetscLineSearchPostCheckFunc func,void *ctx) 235 { 236 PetscFunctionBegin; 237 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_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__ "PetscLineSearchGetPostCheck" 246 /*@C 247 PetscLineSearchGetPostCheck - Gets the post-check function for the line search routine. 248 249 Input Parameters: 250 . linesearch - the PetscLineSearch 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: PetscLineSearchGetPreCheck(), PetscLineSearchSetPostCheck() 262 @*/ 263 PetscErrorCode PetscLineSearchGetPostCheck(PetscLineSearch linesearch, PetscLineSearchPostCheckFunc *func,void **ctx) 264 { 265 PetscFunctionBegin; 266 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_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__ "PetscLineSearchPreCheck" 275 /*@ 276 PetscLineSearchPreCheck - Prepares the line search for being applied. 277 278 Collective on PetscLineSearch 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: PetscLineSearch, Create 289 290 .seealso: PetscLineSearchPostCheck() 291 @*/ 292 PetscErrorCode PetscLineSearchPreCheck(PetscLineSearch 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);CHKERRQ(ierr); 299 } 300 PetscFunctionReturn(0); 301 } 302 303 #undef __FUNCT__ 304 #define __FUNCT__ "PetscLineSearchPostCheck" 305 /*@ 306 PetscLineSearchPostCheck - Prepares the line search for being applied. 307 308 Collective on PetscLineSearch 309 310 Input Parameters: 311 . linesearch - The linesearch instance. 312 313 Output Parameters: 314 + changed_Y - Indicator if the solution has been changed. 315 - changed_W - Indicator if the direction has been changed. 316 317 Level: Intermediate 318 319 .keywords: PetscLineSearch, Create 320 321 .seealso: PetscLineSearchPreCheck() 322 @*/ 323 PetscErrorCode PetscLineSearchPostCheck(PetscLineSearch 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);CHKERRQ(ierr); 331 } 332 PetscFunctionReturn(0); 333 } 334 335 336 #undef __FUNCT__ 337 #define __FUNCT__ "PetscLineSearchPreCheckPicard" 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 PetscLineSearchPreCheckPicard(PetscLineSearch linesearch,Vec X,Vec Y,PetscBool *changed) 371 { 372 PetscErrorCode ierr; 373 PetscReal angle = *(PetscReal*)linesearch->precheckctx; 374 Vec Ylast; 375 PetscScalar dot; 376 PetscInt iter; 377 PetscReal ynorm,ylastnorm,theta,angle_radians; 378 SNES snes; 379 380 PetscFunctionBegin; 381 ierr = PetscLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 382 ierr = PetscObjectQuery((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject*)&Ylast);CHKERRQ(ierr); 383 if (!Ylast) { 384 ierr = VecDuplicate(Y,&Ylast);CHKERRQ(ierr); 385 ierr = PetscObjectCompose((PetscObject)snes,"SNESLineSearchPreCheckPicard_Ylast",(PetscObject)Ylast);CHKERRQ(ierr); 386 ierr = PetscObjectDereference((PetscObject)Ylast);CHKERRQ(ierr); 387 } 388 ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr); 389 if (iter < 2) { 390 ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 391 *changed = PETSC_FALSE; 392 PetscFunctionReturn(0); 393 } 394 395 ierr = VecDot(Y,Ylast,&dot);CHKERRQ(ierr); 396 ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 397 ierr = VecNorm(Ylast,NORM_2,&ylastnorm);CHKERRQ(ierr); 398 /* Compute the angle between the vectors Y and Ylast, clip to keep inside the domain of acos() */ 399 theta = acos((double)PetscClipInterval(PetscAbsScalar(dot) / (ynorm * ylastnorm),-1.0,1.0)); 400 angle_radians = angle * PETSC_PI / 180.; 401 if (PetscAbsReal(theta) < angle_radians || PetscAbsReal(theta - PETSC_PI) < angle_radians) { 402 /* Modify the step Y */ 403 PetscReal alpha,ydiffnorm; 404 ierr = VecAXPY(Ylast,-1.0,Y);CHKERRQ(ierr); 405 ierr = VecNorm(Ylast,NORM_2,&ydiffnorm);CHKERRQ(ierr); 406 alpha = ylastnorm / ydiffnorm; 407 ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 408 ierr = VecScale(Y,alpha);CHKERRQ(ierr); 409 ierr = PetscInfo3(snes,"Angle %G degrees less than threshold %G, corrected step by alpha=%G\n",theta*180./PETSC_PI,angle,alpha);CHKERRQ(ierr); 410 } else { 411 ierr = PetscInfo2(snes,"Angle %G degrees exceeds threshold %G, no correction applied\n",theta*180./PETSC_PI,angle);CHKERRQ(ierr); 412 ierr = VecCopy(Y,Ylast);CHKERRQ(ierr); 413 *changed = PETSC_FALSE; 414 } 415 PetscFunctionReturn(0); 416 } 417 418 #undef __FUNCT__ 419 #define __FUNCT__ "PetscLineSearchApply" 420 /*@ 421 PetscLineSearchApply - Computes the line-search update 422 423 Collective on PetscLineSearch 424 425 Input Parameters: 426 + linesearch - The linesearch instance. 427 . X - The current solution. 428 . F - The current function. 429 . fnorm - The current norm. 430 . Y - The search direction. 431 432 Output Parameters: 433 + X - The new solution. 434 . F - The new function. 435 - fnorm - The new function norm. 436 437 Level: Intermediate 438 439 .keywords: PetscLineSearch, Create 440 441 .seealso: PetscLineSearchCreate(), PetscLineSearchPreCheck(), PetscLineSearchPostCheck() 442 @*/ 443 PetscErrorCode PetscLineSearchApply(PetscLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) { 444 PetscErrorCode ierr; 445 PetscFunctionBegin; 446 447 /* check the pointers */ 448 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 449 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 450 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 451 PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 452 453 linesearch->success = PETSC_TRUE; 454 455 linesearch->vec_sol = X; 456 linesearch->vec_update = Y; 457 linesearch->vec_func = F; 458 459 ierr = PetscLineSearchSetUp(linesearch);CHKERRQ(ierr); 460 461 if (!linesearch->keeplambda) 462 linesearch->lambda = linesearch->damping; /* set the initial guess to lambda */ 463 464 if (fnorm) { 465 linesearch->fnorm = *fnorm; 466 } else { 467 ierr = VecNorm(F, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 468 } 469 470 ierr = PetscLogEventBegin(PetscLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 471 472 ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr); 473 474 ierr = PetscLogEventEnd(PetscLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 475 476 if (fnorm) 477 *fnorm = linesearch->fnorm; 478 PetscFunctionReturn(0); 479 } 480 481 #undef __FUNCT__ 482 #define __FUNCT__ "PetscLineSearchDestroy" 483 /*@ 484 PetscLineSearchDestroy - Destroys the line search instance. 485 486 Collective on PetscLineSearch 487 488 Input Parameters: 489 . linesearch - The linesearch instance. 490 491 Level: Intermediate 492 493 .keywords: PetscLineSearch, Create 494 495 .seealso: PetscLineSearchCreate(), PetscLineSearchReset() 496 @*/ 497 PetscErrorCode PetscLineSearchDestroy(PetscLineSearch * linesearch) { 498 PetscErrorCode ierr; 499 PetscFunctionBegin; 500 if (!*linesearch) PetscFunctionReturn(0); 501 PetscValidHeaderSpecific((*linesearch),PETSCLINESEARCH_CLASSID,1); 502 if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);} 503 ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr); 504 ierr = PetscLineSearchReset(*linesearch); 505 if ((*linesearch)->ops->destroy) { 506 (*linesearch)->ops->destroy(*linesearch); 507 } 508 ierr = PetscViewerDestroy(&(*linesearch)->monitor);CHKERRQ(ierr); 509 ierr = PetscHeaderDestroy(linesearch);CHKERRQ(ierr); 510 PetscFunctionReturn(0); 511 } 512 513 #undef __FUNCT__ 514 #define __FUNCT__ "PetscLineSearchSetMonitor" 515 /*@ 516 PetscLineSearchSetMonitor - Turns on/off printing useful things about the line search. 517 518 Input Parameters: 519 + snes - nonlinear context obtained from SNESCreate() 520 - flg - PETSC_TRUE to monitor the line search 521 522 Logically Collective on SNES 523 524 Options Database: 525 . -linesearch_monitor - enables the monitor. 526 527 Level: intermediate 528 529 530 .seealso: PetscLineSearchGetMonitor() 531 @*/ 532 PetscErrorCode PetscLineSearchSetMonitor(PetscLineSearch linesearch, PetscBool flg) 533 { 534 535 PetscErrorCode ierr; 536 PetscFunctionBegin; 537 if (flg && !linesearch->monitor) { 538 ierr = PetscViewerASCIIOpen(((PetscObject)linesearch)->comm,"stdout",&linesearch->monitor);CHKERRQ(ierr); 539 } else if (!flg && linesearch->monitor) { 540 ierr = PetscViewerDestroy(&linesearch->monitor);CHKERRQ(ierr); 541 } 542 PetscFunctionReturn(0); 543 } 544 545 #undef __FUNCT__ 546 #define __FUNCT__ "PetscLineSearchGetMonitor" 547 /*@ 548 PetscLineSearchGetMonitor - Gets the monitor instance for the line search 549 550 Input Parameters: 551 . linesearch - linesearch context. 552 553 Input Parameters: 554 . monitor - monitor context. 555 556 Logically Collective on SNES 557 558 559 Options Database Keys: 560 . -linesearch_monitor - enables the monitor. 561 562 Level: intermediate 563 564 565 .seealso: PetscLineSearchSetMonitor() 566 @*/ 567 PetscErrorCode PetscLineSearchGetMonitor(PetscLineSearch linesearch, PetscViewer *monitor) 568 { 569 570 PetscFunctionBegin; 571 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 572 if (monitor) { 573 PetscValidPointer(monitor, 2); 574 *monitor = linesearch->monitor; 575 } 576 PetscFunctionReturn(0); 577 } 578 579 #undef __FUNCT__ 580 #define __FUNCT__ "PetscLineSearchSetFromOptions" 581 /*@ 582 PetscLineSearchSetFromOptions - Sets options for the line search 583 584 Input Parameters: 585 . linesearch - linesearch context. 586 587 Options Database Keys: 588 + -linesearch_type - The Line search method 589 . -linesearch_monitor - Print progress of line searches 590 . -linesearch_damping - The linesearch damping parameter. 591 . -linesearch_norms - Turn on/off the linesearch norms 592 . -linesearch_keeplambda - Keep the previous search length as the initial guess. 593 - -linesearch_max_it - The number of iterations for iterative line searches. 594 595 Logically Collective on PetscLineSearch 596 597 Level: intermediate 598 599 600 .seealso: PetscLineSearchCreate() 601 @*/ 602 PetscErrorCode PetscLineSearchSetFromOptions(PetscLineSearch linesearch) { 603 PetscErrorCode ierr; 604 const char *deft = PETSCLINESEARCHBASIC; 605 char type[256]; 606 PetscBool flg, set; 607 PetscFunctionBegin; 608 if (!PetscLineSearchRegisterAllCalled) {ierr = PetscLineSearchRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 609 610 ierr = PetscObjectOptionsBegin((PetscObject)linesearch);CHKERRQ(ierr); 611 if (((PetscObject)linesearch)->type_name) { 612 deft = ((PetscObject)linesearch)->type_name; 613 } 614 ierr = PetscOptionsList("-linesearch_type","Line-search method","PetscLineSearchSetType",PetscLineSearchList,deft,type,256,&flg);CHKERRQ(ierr); 615 if (flg) { 616 ierr = PetscLineSearchSetType(linesearch,type);CHKERRQ(ierr); 617 } else if (!((PetscObject)linesearch)->type_name) { 618 ierr = PetscLineSearchSetType(linesearch,deft);CHKERRQ(ierr); 619 } 620 if (linesearch->ops->setfromoptions) { 621 (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr); 622 } 623 624 ierr = PetscOptionsBool("-linesearch_monitor","Print progress of line searches","SNESPetscLineSearchSetMonitor", 625 linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 626 if (set) {ierr = PetscLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);} 627 628 ierr = PetscOptionsReal("-linesearch_damping","Line search damping and initial step guess","PetscLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr); 629 ierr = PetscOptionsReal("-linesearch_rtol","Tolerance for iterative line search","PetscLineSearchSetRTolerance",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr); 630 ierr = PetscOptionsBool("-linesearch_norms","Compute final norms in line search","PetscLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr); 631 ierr = PetscOptionsBool("-linesearch_keeplambda","Use previous lambda as damping","PetscLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr); 632 ierr = PetscOptionsInt("-linesearch_max_it","Maximum iterations for iterative line searches","",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr); 633 ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr); 634 ierr = PetscOptionsEnd();CHKERRQ(ierr); 635 PetscFunctionReturn(0); 636 } 637 638 #undef __FUNCT__ 639 #define __FUNCT__ "PetscLineSearchView" 640 /*@ 641 PetscLineSearchView - Views useful information for the line search. 642 643 Input Parameters: 644 . linesearch - linesearch context. 645 646 Logically Collective on PetscLineSearch 647 648 Level: intermediate 649 650 651 .seealso: PetscLineSearchCreate() 652 @*/ 653 PetscErrorCode PetscLineSearchView(PetscLineSearch linesearch) { 654 PetscFunctionBegin; 655 656 PetscFunctionReturn(0); 657 } 658 659 #undef __FUNCT__ 660 #define __FUNCT__ "PetscLineSearchSetType" 661 /*@C 662 PetscLineSearchSetType - Sets the linesearch type 663 664 Input Parameters: 665 + linesearch - linesearch context. 666 - type - The type of line search to be used 667 668 Logically Collective on PetscLineSearch 669 670 Level: intermediate 671 672 673 .seealso: PetscLineSearchCreate() 674 @*/ 675 PetscErrorCode PetscLineSearchSetType(PetscLineSearch linesearch, const PetscLineSearchType type) 676 { 677 678 PetscErrorCode ierr,(*r)(PetscLineSearch); 679 PetscBool match; 680 681 PetscFunctionBegin; 682 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 683 PetscValidCharPointer(type,2); 684 685 ierr = PetscTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr); 686 if (match) PetscFunctionReturn(0); 687 688 ierr = PetscFListFind(PetscLineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 689 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type); 690 /* Destroy the previous private linesearch context */ 691 if (linesearch->ops->destroy) { 692 ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr); 693 linesearch->ops->destroy = PETSC_NULL; 694 } 695 /* Reinitialize function pointers in PetscLineSearchOps structure */ 696 linesearch->ops->apply = 0; 697 linesearch->ops->view = 0; 698 linesearch->ops->setfromoptions = 0; 699 linesearch->ops->destroy = 0; 700 701 ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr); 702 ierr = (*r)(linesearch);CHKERRQ(ierr); 703 #if defined(PETSC_HAVE_AMS) 704 if (PetscAMSPublishAll) { 705 ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr); 706 } 707 #endif 708 PetscFunctionReturn(0); 709 } 710 711 #undef __FUNCT__ 712 #define __FUNCT__ "PetscLineSearchSetSNES" 713 /*@ 714 PetscLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation 715 716 Input Parameters: 717 + linesearch - linesearch context. 718 - snes - The snes instance 719 720 Level: intermediate 721 722 723 .seealso: PetscLineSearchGetSNES(), PetscLineSearchSetVecs() 724 @*/ 725 PetscErrorCode PetscLineSearchSetSNES(PetscLineSearch linesearch, SNES snes){ 726 PetscFunctionBegin; 727 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 728 PetscValidHeaderSpecific(snes,SNES_CLASSID,2); 729 linesearch->snes = snes; 730 PetscFunctionReturn(0); 731 } 732 733 #undef __FUNCT__ 734 #define __FUNCT__ "PetscLineSearchGetSNES" 735 /*@ 736 PetscLineSearchGetSNES - Gets the SNES for the linesearch for function evaluation 737 738 Input Parameters: 739 . linesearch - linesearch context. 740 741 Output Parameters: 742 . snes - The snes instance 743 744 Level: intermediate 745 746 .seealso: PetscLineSearchGetSNES(), PetscLineSearchSetVecs() 747 @*/ 748 PetscErrorCode PetscLineSearchGetSNES(PetscLineSearch linesearch, SNES *snes){ 749 PetscFunctionBegin; 750 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 751 PetscValidPointer(snes, 2); 752 *snes = linesearch->snes; 753 PetscFunctionReturn(0); 754 } 755 756 #undef __FUNCT__ 757 #define __FUNCT__ "PetscLineSearchGetLambda" 758 /*@ 759 PetscLineSearchGetLambda - Gets the last linesearch steplength discovered. 760 761 Input Parameters: 762 . linesearch - linesearch context. 763 764 Output Parameters: 765 . lambda - The last steplength. 766 767 Level: intermediate 768 769 .seealso: PetscLineSearchGetSNES(), PetscLineSearchSetVecs() 770 @*/ 771 PetscErrorCode PetscLineSearchGetLambda(PetscLineSearch linesearch,PetscReal *lambda) 772 { 773 PetscFunctionBegin; 774 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 775 PetscValidPointer(lambda, 2); 776 *lambda = linesearch->lambda; 777 PetscFunctionReturn(0); 778 } 779 780 #undef __FUNCT__ 781 #define __FUNCT__ "PetscLineSearchSetLambda" 782 /*@ 783 PetscLineSearchSetLambda - Sets the linesearch steplength. 784 785 Input Parameters: 786 + linesearch - linesearch context. 787 - lambda - The last steplength. 788 789 Level: intermediate 790 791 .seealso: PetscLineSearchGetLambda() 792 @*/ 793 PetscErrorCode PetscLineSearchSetLambda(PetscLineSearch linesearch, PetscReal lambda) 794 { 795 PetscFunctionBegin; 796 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 797 linesearch->lambda = lambda; 798 PetscFunctionReturn(0); 799 } 800 801 #undef __FUNCT__ 802 #define __FUNCT__ "PetscLineSearchGetTolerances" 803 /*@ 804 PetscLineSearchGetTolerances - Gets the tolerances for the method 805 806 Input Parameters: 807 . linesearch - linesearch context. 808 809 Output Parameters: 810 + steptol - The minimum steplength 811 . rtol - The relative tolerance for iterative line searches 812 . atol - The absolute tolerance for iterative line searches 813 . ltol - The change in lambda tolerance for iterative line searches 814 - max_it - The maximum number of iterations of the line search 815 816 817 Level: advanced 818 819 .seealso: PetscLineSearchSetTolerances() 820 @*/ 821 PetscErrorCode PetscLineSearchGetTolerances(PetscLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 822 { 823 PetscFunctionBegin; 824 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 825 if (steptol) { 826 PetscValidPointer(steptol, 2); 827 *steptol = linesearch->steptol; 828 } 829 if (maxstep) { 830 PetscValidPointer(maxstep, 3); 831 *maxstep = linesearch->maxstep; 832 } 833 if (rtol) { 834 PetscValidPointer(rtol, 4); 835 *rtol = linesearch->rtol; 836 } 837 if (atol) { 838 PetscValidPointer(atol, 5); 839 *atol = linesearch->atol; 840 } 841 if (ltol) { 842 PetscValidPointer(ltol, 6); 843 *ltol = linesearch->ltol; 844 } 845 if (max_its) { 846 PetscValidPointer(max_its, 7); 847 *max_its = linesearch->max_its; 848 } 849 PetscFunctionReturn(0); 850 } 851 852 #undef __FUNCT__ 853 #define __FUNCT__ "PetscLineSearchSetTolerances" 854 /*@ 855 PetscLineSearchSetTolerances - Sets the tolerances for the method 856 857 Input Parameters: 858 + linesearch - linesearch context. 859 . steptol - The minimum steplength 860 . rtol - The relative tolerance for iterative line searches 861 . atol - The absolute tolerance for iterative line searches 862 . ltol - The change in lambda tolerance for iterative line searches 863 - max_it - The maximum number of iterations of the line search 864 865 866 Level: advanced 867 868 .seealso: PetscLineSearchGetTolerances() 869 @*/ 870 PetscErrorCode PetscLineSearchSetTolerances(PetscLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 871 { 872 PetscFunctionBegin; 873 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 874 linesearch->steptol = steptol; 875 linesearch->maxstep = maxstep; 876 linesearch->rtol = rtol; 877 linesearch->atol = atol; 878 linesearch->ltol = ltol; 879 linesearch->max_its = max_its; 880 PetscFunctionReturn(0); 881 } 882 883 884 #undef __FUNCT__ 885 #define __FUNCT__ "PetscLineSearchGetDamping" 886 /*@ 887 PetscLineSearchGetDamping - Gets the line search damping parameter. 888 889 Input Parameters: 890 . linesearch - linesearch context. 891 892 Output Parameters: 893 . damping - The damping parameter. 894 895 Level: intermediate 896 897 .seealso: PetscLineSearchGetStepTolerance() 898 @*/ 899 900 PetscErrorCode PetscLineSearchGetDamping(PetscLineSearch linesearch,PetscReal *damping) 901 { 902 PetscFunctionBegin; 903 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 904 PetscValidPointer(damping, 2); 905 *damping = linesearch->damping; 906 PetscFunctionReturn(0); 907 } 908 909 #undef __FUNCT__ 910 #define __FUNCT__ "PetscLineSearchSetDamping" 911 /*@ 912 PetscLineSearchSetDamping - Sets the line search damping paramter. 913 914 Input Parameters: 915 . linesearch - linesearch context. 916 . damping - The damping parameter. 917 918 Level: intermediate 919 920 .seealso: PetscLineSearchGetDamping() 921 @*/ 922 PetscErrorCode PetscLineSearchSetDamping(PetscLineSearch linesearch,PetscReal damping) 923 { 924 PetscFunctionBegin; 925 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 926 linesearch->damping = damping; 927 PetscFunctionReturn(0); 928 } 929 930 #undef __FUNCT__ 931 #define __FUNCT__ "PetscLineSearchGetNorms" 932 /*@ 933 PetscLineSearchGetNorms - Gets the norms for for X, Y, and F. 934 935 Input Parameters: 936 . linesearch - linesearch context. 937 938 Output Parameters: 939 + xnorm - The norm of the current solution 940 . fnorm - The norm of the current function 941 - ynorm - The norm of the current update 942 943 Level: intermediate 944 945 .seealso: PetscLineSearchSetNorms() PetscLineSearchGetVecs() 946 @*/ 947 PetscErrorCode PetscLineSearchGetNorms(PetscLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm) 948 { 949 PetscFunctionBegin; 950 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 951 if (xnorm) { 952 *xnorm = linesearch->xnorm; 953 } 954 if (fnorm) { 955 *fnorm = linesearch->fnorm; 956 } 957 if (ynorm) { 958 *ynorm = linesearch->ynorm; 959 } 960 PetscFunctionReturn(0); 961 } 962 963 #undef __FUNCT__ 964 #define __FUNCT__ "PetscLineSearchSetNorms" 965 /*@ 966 PetscLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 967 968 Input Parameters: 969 + linesearch - linesearch context. 970 . xnorm - The norm of the current solution 971 . fnorm - The norm of the current function 972 - ynorm - The norm of the current update 973 974 Level: intermediate 975 976 .seealso: PetscLineSearchGetNorms(), PetscLineSearchSetVecs() 977 @*/ 978 PetscErrorCode PetscLineSearchSetNorms(PetscLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 979 { 980 PetscFunctionBegin; 981 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 982 linesearch->xnorm = xnorm; 983 linesearch->fnorm = fnorm; 984 linesearch->ynorm = ynorm; 985 PetscFunctionReturn(0); 986 } 987 988 #undef __FUNCT__ 989 #define __FUNCT__ "PetscLineSearchComputeNorms" 990 /*@ 991 PetscLineSearchComputeNorms - Computes the norms of X, F, and Y. 992 993 Input Parameters: 994 . linesearch - linesearch context. 995 996 Options Database Keys: 997 . -linesearch_norms - turn norm computation on or off. 998 999 Level: intermediate 1000 1001 .seealso: PetscLineSearchGetNorms, PetscLineSearchSetNorms() 1002 @*/ 1003 PetscErrorCode PetscLineSearchComputeNorms(PetscLineSearch linesearch) 1004 { 1005 PetscErrorCode ierr; 1006 SNES snes; 1007 PetscFunctionBegin; 1008 if (linesearch->norms) { 1009 if (linesearch->ops->vinorm) { 1010 ierr = PetscLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 1011 ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1012 ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1013 ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr); 1014 } else { 1015 ierr = VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 1016 ierr = VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1017 ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1018 ierr = VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 1019 ierr = VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1020 ierr = VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1021 } 1022 } 1023 PetscFunctionReturn(0); 1024 } 1025 1026 #undef __FUNCT__ 1027 #define __FUNCT__ "PetscLineSearchGetVecs" 1028 /*@ 1029 PetscLineSearchGetVecs - Gets the vectors from the PetscLineSearch context 1030 1031 Input Parameters: 1032 . linesearch - linesearch context. 1033 1034 Output Parameters: 1035 + X - The old solution 1036 . F - The old function 1037 . Y - The search direction 1038 . W - The new solution 1039 - G - The new function 1040 1041 Level: intermediate 1042 1043 .seealso: PetscLineSearchGetNorms(), PetscLineSearchSetVecs() 1044 @*/ 1045 PetscErrorCode PetscLineSearchGetVecs(PetscLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) { 1046 PetscFunctionBegin; 1047 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 1048 if (X) { 1049 PetscValidPointer(X, 2); 1050 *X = linesearch->vec_sol; 1051 } 1052 if (F) { 1053 PetscValidPointer(F, 3); 1054 *F = linesearch->vec_func; 1055 } 1056 if (Y) { 1057 PetscValidPointer(Y, 4); 1058 *Y = linesearch->vec_update; 1059 } 1060 if (W) { 1061 PetscValidPointer(W, 5); 1062 *W = linesearch->vec_sol_new; 1063 } 1064 if (G) { 1065 PetscValidPointer(G, 6); 1066 *G = linesearch->vec_func_new; 1067 } 1068 1069 PetscFunctionReturn(0); 1070 } 1071 1072 #undef __FUNCT__ 1073 #define __FUNCT__ "PetscLineSearchSetVecs" 1074 /*@ 1075 PetscLineSearchSetVecs - Sets the vectors on the PetscLineSearch context 1076 1077 Input Parameters: 1078 + linesearch - linesearch context. 1079 . X - The old solution 1080 . F - The old function 1081 . Y - The search direction 1082 . W - The new solution 1083 - G - The new function 1084 1085 Level: intermediate 1086 1087 .seealso: PetscLineSearchSetNorms(), PetscLineSearchGetVecs() 1088 @*/ 1089 PetscErrorCode PetscLineSearchSetVecs(PetscLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) { 1090 PetscFunctionBegin; 1091 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 1092 if (X) { 1093 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1094 linesearch->vec_sol = X; 1095 } 1096 if (F) { 1097 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 1098 linesearch->vec_func = F; 1099 } 1100 if (Y) { 1101 PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 1102 linesearch->vec_update = Y; 1103 } 1104 if (W) { 1105 PetscValidHeaderSpecific(W,VEC_CLASSID,5); 1106 linesearch->vec_sol_new = W; 1107 } 1108 if (G) { 1109 PetscValidHeaderSpecific(G,VEC_CLASSID,6); 1110 linesearch->vec_func_new = G; 1111 } 1112 1113 PetscFunctionReturn(0); 1114 } 1115 1116 #undef __FUNCT__ 1117 #define __FUNCT__ "PetscLineSearchAppendOptionsPrefix" 1118 /*@C 1119 PetscLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1120 SNES options in the database. 1121 1122 Logically Collective on SNES 1123 1124 Input Parameters: 1125 + snes - the SNES context 1126 - prefix - the prefix to prepend to all option names 1127 1128 Notes: 1129 A hyphen (-) must NOT be given at the beginning of the prefix name. 1130 The first character of all runtime options is AUTOMATICALLY the hyphen. 1131 1132 Level: advanced 1133 1134 .keywords: PetscLineSearch, append, options, prefix, database 1135 1136 .seealso: SNESGetOptionsPrefix() 1137 @*/ 1138 PetscErrorCode PetscLineSearchAppendOptionsPrefix(PetscLineSearch linesearch,const char prefix[]) 1139 { 1140 PetscErrorCode ierr; 1141 1142 PetscFunctionBegin; 1143 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 1144 ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1145 PetscFunctionReturn(0); 1146 } 1147 1148 #undef __FUNCT__ 1149 #define __FUNCT__ "PetscLineSearchGetOptionsPrefix" 1150 /*@C 1151 PetscLineSearchGetOptionsPrefix - Sets the prefix used for searching for all 1152 PetscLineSearch options in the database. 1153 1154 Not Collective 1155 1156 Input Parameter: 1157 . snes - the SNES context 1158 1159 Output Parameter: 1160 . prefix - pointer to the prefix string used 1161 1162 Notes: On the fortran side, the user should pass in a string 'prefix' of 1163 sufficient length to hold the prefix. 1164 1165 Level: advanced 1166 1167 .keywords: PetscLineSearch, get, options, prefix, database 1168 1169 .seealso: SNESAppendOptionsPrefix() 1170 @*/ 1171 PetscErrorCode PetscLineSearchGetOptionsPrefix(PetscLineSearch linesearch,const char *prefix[]) 1172 { 1173 PetscErrorCode ierr; 1174 1175 PetscFunctionBegin; 1176 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 1177 ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1178 PetscFunctionReturn(0); 1179 } 1180 1181 #undef __FUNCT__ 1182 #define __FUNCT__ "PetscLineSearchGetWork" 1183 /*@ 1184 PetscLineSearchGetWork - Gets work vectors for the line search. 1185 1186 Input Parameter: 1187 + linesearch - the PetscLineSearch context 1188 - nwork - the number of work vectors 1189 1190 Level: developer 1191 1192 .keywords: PetscLineSearch, work, vector 1193 1194 .seealso: SNESDefaultGetWork() 1195 @*/ 1196 PetscErrorCode PetscLineSearchGetWork(PetscLineSearch linesearch, PetscInt nwork) 1197 { 1198 PetscErrorCode ierr; 1199 PetscFunctionBegin; 1200 if (linesearch->vec_sol) { 1201 ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr); 1202 } else { 1203 SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 1204 } 1205 PetscFunctionReturn(0); 1206 } 1207 1208 1209 #undef __FUNCT__ 1210 #define __FUNCT__ "PetscLineSearchGetSuccess" 1211 /*@ 1212 PetscLineSearchGetSuccess - Gets the success/failure status of the last line search application 1213 1214 Input Parameters: 1215 . linesearch - linesearch context. 1216 1217 Output Parameters: 1218 . success - The success or failure status. 1219 1220 Level: intermediate 1221 1222 .seealso: PetscLineSearchSetSuccess() 1223 @*/ 1224 PetscErrorCode PetscLineSearchGetSuccess(PetscLineSearch linesearch, PetscBool *success) 1225 { 1226 PetscFunctionBegin; 1227 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 1228 PetscValidPointer(success, 2); 1229 if (success) { 1230 *success = linesearch->success; 1231 } 1232 PetscFunctionReturn(0); 1233 } 1234 1235 #undef __FUNCT__ 1236 #define __FUNCT__ "PetscLineSearchSetSuccess" 1237 /*@ 1238 PetscLineSearchSetSuccess - Sets the success/failure status of the last line search application 1239 1240 Input Parameters: 1241 + linesearch - linesearch context. 1242 - success - The success or failure status. 1243 1244 Level: intermediate 1245 1246 .seealso: PetscLineSearchGetSuccess() 1247 @*/ 1248 PetscErrorCode PetscLineSearchSetSuccess(PetscLineSearch linesearch, PetscBool success) 1249 { 1250 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 1251 PetscFunctionBegin; 1252 linesearch->success = success; 1253 PetscFunctionReturn(0); 1254 } 1255 1256 #undef __FUNCT__ 1257 #define __FUNCT__ "PetscLineSearchSetVIFunctions" 1258 /*@C 1259 PetscLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 1260 1261 Input Parameters: 1262 + snes - nonlinear context obtained from SNESCreate() 1263 . projectfunc - function for projecting the function to the bounds 1264 - normfunc - function for computing the norm of an active set 1265 1266 Logically Collective on SNES 1267 1268 Calling sequence of projectfunc: 1269 .vb 1270 projectfunc (SNES snes, Vec X) 1271 .ve 1272 1273 Input parameters for projectfunc: 1274 + snes - nonlinear context 1275 - X - current solution 1276 1277 Output parameters for func: 1278 . X - Projected solution 1279 1280 Calling sequence of normfunc: 1281 .vb 1282 projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm) 1283 .ve 1284 1285 Input parameters for projectfunc: 1286 + snes - nonlinear context 1287 . X - current solution 1288 - F - current residual 1289 1290 Output parameters for func: 1291 . fnorm - VI-specific norm of the function 1292 1293 1294 Level: developer 1295 1296 .keywords: SNES, line search, VI, nonlinear, set, line search 1297 1298 .seealso: PetscLineSearchGetVIFunctions(), PetscLineSearchSetPostCheck(), PetscLineSearchSetPreCheck() 1299 @*/ 1300 extern PetscErrorCode PetscLineSearchSetVIFunctions(PetscLineSearch linesearch, PetscLineSearchVIProjectFunc projectfunc, PetscLineSearchVINormFunc normfunc) 1301 { 1302 PetscFunctionBegin; 1303 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 1304 if (projectfunc) linesearch->ops->viproject = projectfunc; 1305 if (normfunc) linesearch->ops->vinorm = normfunc; 1306 PetscFunctionReturn(0); 1307 } 1308 1309 #undef __FUNCT__ 1310 #define __FUNCT__ "PetscLineSearchGetVIFunctions" 1311 /*@C 1312 PetscLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 1313 1314 Input Parameters: 1315 . snes - nonlinear context obtained from SNESCreate() 1316 1317 Output Parameters: 1318 + projectfunc - function for projecting the function to the bounds 1319 - normfunc - function for computing the norm of an active set 1320 1321 Logically Collective on SNES 1322 1323 Level: developer 1324 1325 .keywords: SNES, line search, VI, nonlinear, get, line search 1326 1327 .seealso: PetscLineSearchSetVIFunctions(), PetscLineSearchGetPostCheck(), PetscLineSearchGetPreCheck() 1328 @*/ 1329 extern PetscErrorCode PetscLineSearchGetVIFunctions(PetscLineSearch linesearch, PetscLineSearchVIProjectFunc *projectfunc, PetscLineSearchVINormFunc *normfunc) 1330 { 1331 PetscFunctionBegin; 1332 if (projectfunc) *projectfunc = linesearch->ops->viproject; 1333 if (normfunc) *normfunc = linesearch->ops->vinorm; 1334 PetscFunctionReturn(0); 1335 } 1336 1337 #undef __FUNCT__ 1338 #define __FUNCT__ "PetscLineSearchRegister" 1339 /*@C 1340 PetscLineSearchRegister - See PetscLineSearchRegisterDynamic() 1341 1342 Level: advanced 1343 @*/ 1344 PetscErrorCode PetscLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(PetscLineSearch)) 1345 { 1346 char fullname[PETSC_MAX_PATH_LEN]; 1347 PetscErrorCode ierr; 1348 1349 PetscFunctionBegin; 1350 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 1351 ierr = PetscFListAdd(&PetscLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 1352 PetscFunctionReturn(0); 1353 } 1354