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, linesearch->precheckctx);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, linesearch->postcheckctx);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,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 = PetscLineSearchGetSNES(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__ "PetscLineSearchApply" 421 /*@ 422 PetscLineSearchApply - Computes the line-search update 423 424 Collective on PetscLineSearch 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: PetscLineSearch, Create 441 442 .seealso: PetscLineSearchCreate(), PetscLineSearchPreCheck(), PetscLineSearchPostCheck() 443 @*/ 444 PetscErrorCode PetscLineSearchApply(PetscLineSearch linesearch, Vec X, Vec F, PetscReal * fnorm, Vec Y) { 445 PetscErrorCode ierr; 446 PetscFunctionBegin; 447 448 /* check the pointers */ 449 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_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 = PetscLineSearchSetUp(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(PetscLineSearch_Apply,linesearch,X,F,Y);CHKERRQ(ierr); 472 473 ierr = (*linesearch->ops->apply)(linesearch);CHKERRQ(ierr); 474 475 ierr = PetscLogEventEnd(PetscLineSearch_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__ "PetscLineSearchDestroy" 484 /*@ 485 PetscLineSearchDestroy - Destroys the line search instance. 486 487 Collective on PetscLineSearch 488 489 Input Parameters: 490 . linesearch - The linesearch instance. 491 492 Level: Intermediate 493 494 .keywords: PetscLineSearch, Create 495 496 .seealso: PetscLineSearchCreate(), PetscLineSearchReset() 497 @*/ 498 PetscErrorCode PetscLineSearchDestroy(PetscLineSearch * linesearch) { 499 PetscErrorCode ierr; 500 PetscFunctionBegin; 501 if (!*linesearch) PetscFunctionReturn(0); 502 PetscValidHeaderSpecific((*linesearch),PETSCLINESEARCH_CLASSID,1); 503 if (--((PetscObject)(*linesearch))->refct > 0) {*linesearch = 0; PetscFunctionReturn(0);} 504 ierr = PetscObjectDepublish((*linesearch));CHKERRQ(ierr); 505 ierr = PetscLineSearchReset(*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__ "PetscLineSearchSetMonitor" 516 /*@ 517 PetscLineSearchSetMonitor - 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: PetscLineSearchGetMonitor() 532 @*/ 533 PetscErrorCode PetscLineSearchSetMonitor(PetscLineSearch 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__ "PetscLineSearchGetMonitor" 548 /*@ 549 PetscLineSearchGetMonitor - 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: PetscLineSearchSetMonitor() 567 @*/ 568 PetscErrorCode PetscLineSearchGetMonitor(PetscLineSearch linesearch, PetscViewer *monitor) 569 { 570 571 PetscFunctionBegin; 572 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_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__ "PetscLineSearchSetFromOptions" 582 /*@ 583 PetscLineSearchSetFromOptions - 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 PetscLineSearch 597 598 Level: intermediate 599 600 601 .seealso: PetscLineSearchCreate() 602 @*/ 603 PetscErrorCode PetscLineSearchSetFromOptions(PetscLineSearch linesearch) { 604 PetscErrorCode ierr; 605 const char *deft = PETSCLINESEARCHBASIC; 606 char type[256]; 607 PetscBool flg, set; 608 PetscFunctionBegin; 609 if (!PetscLineSearchRegisterAllCalled) {ierr = PetscLineSearchRegisterAll(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("-linesearch_type","Line-search method","PetscLineSearchSetType",PetscLineSearchList,deft,type,256,&flg);CHKERRQ(ierr); 616 if (flg) { 617 ierr = PetscLineSearchSetType(linesearch,type);CHKERRQ(ierr); 618 } else if (!((PetscObject)linesearch)->type_name) { 619 ierr = PetscLineSearchSetType(linesearch,deft);CHKERRQ(ierr); 620 } 621 if (linesearch->ops->setfromoptions) { 622 (*linesearch->ops->setfromoptions)(linesearch);CHKERRQ(ierr); 623 } 624 625 ierr = PetscOptionsBool("-linesearch_monitor","Print progress of line searches","SNESPetscLineSearchSetMonitor", 626 linesearch->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 627 if (set) {ierr = PetscLineSearchSetMonitor(linesearch,flg);CHKERRQ(ierr);} 628 629 ierr = PetscOptionsReal("-linesearch_damping","Line search damping and initial step guess","PetscLineSearchSetDamping",linesearch->damping,&linesearch->damping,0);CHKERRQ(ierr); 630 ierr = PetscOptionsReal("-linesearch_rtol","Tolerance for iterative line search","PetscLineSearchSetRTolerance",linesearch->rtol,&linesearch->rtol,0);CHKERRQ(ierr); 631 ierr = PetscOptionsBool("-linesearch_norms","Compute final norms in line search","PetscLineSearchSetComputeNorms",linesearch->norms,&linesearch->norms,0);CHKERRQ(ierr); 632 ierr = PetscOptionsBool("-linesearch_keeplambda","Use previous lambda as damping","PetscLineSearchSetKeepLambda",linesearch->keeplambda,&linesearch->keeplambda,0);CHKERRQ(ierr); 633 ierr = PetscOptionsInt("-linesearch_max_it","Maximum iterations for iterative line searches","",linesearch->max_its,&linesearch->max_its,0);CHKERRQ(ierr); 634 ierr = PetscObjectProcessOptionsHandlers((PetscObject)linesearch);CHKERRQ(ierr); 635 ierr = PetscOptionsEnd();CHKERRQ(ierr); 636 PetscFunctionReturn(0); 637 } 638 639 #undef __FUNCT__ 640 #define __FUNCT__ "PetscLineSearchView" 641 /*@ 642 PetscLineSearchView - Views useful information for the line search. 643 644 Input Parameters: 645 . linesearch - linesearch context. 646 647 Logically Collective on PetscLineSearch 648 649 Level: intermediate 650 651 652 .seealso: PetscLineSearchCreate() 653 @*/ 654 PetscErrorCode PetscLineSearchView(PetscLineSearch linesearch) { 655 PetscFunctionBegin; 656 657 PetscFunctionReturn(0); 658 } 659 660 #undef __FUNCT__ 661 #define __FUNCT__ "PetscLineSearchSetType" 662 /*@C 663 PetscLineSearchSetType - Sets the linesearch type 664 665 Input Parameters: 666 + linesearch - linesearch context. 667 - type - The type of line search to be used 668 669 Logically Collective on PetscLineSearch 670 671 Level: intermediate 672 673 674 .seealso: PetscLineSearchCreate() 675 @*/ 676 PetscErrorCode PetscLineSearchSetType(PetscLineSearch linesearch, const PetscLineSearchType type) 677 { 678 679 PetscErrorCode ierr,(*r)(PetscLineSearch); 680 PetscBool match; 681 682 PetscFunctionBegin; 683 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 684 PetscValidCharPointer(type,2); 685 686 ierr = PetscTypeCompare((PetscObject)linesearch,type,&match);CHKERRQ(ierr); 687 if (match) PetscFunctionReturn(0); 688 689 ierr = PetscFListFind(PetscLineSearchList,((PetscObject)linesearch)->comm,type,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 690 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Line Search type %s",type); 691 /* Destroy the previous private linesearch context */ 692 if (linesearch->ops->destroy) { 693 ierr = (*(linesearch)->ops->destroy)(linesearch);CHKERRQ(ierr); 694 linesearch->ops->destroy = PETSC_NULL; 695 } 696 /* Reinitialize function pointers in PetscLineSearchOps structure */ 697 linesearch->ops->apply = 0; 698 linesearch->ops->view = 0; 699 linesearch->ops->setfromoptions = 0; 700 linesearch->ops->destroy = 0; 701 702 ierr = PetscObjectChangeTypeName((PetscObject)linesearch,type);CHKERRQ(ierr); 703 ierr = (*r)(linesearch);CHKERRQ(ierr); 704 #if defined(PETSC_HAVE_AMS) 705 if (PetscAMSPublishAll) { 706 ierr = PetscObjectAMSPublish((PetscObject)linesearch);CHKERRQ(ierr); 707 } 708 #endif 709 PetscFunctionReturn(0); 710 } 711 712 #undef __FUNCT__ 713 #define __FUNCT__ "PetscLineSearchSetSNES" 714 /*@ 715 PetscLineSearchSetSNES - Sets the SNES for the linesearch for function evaluation 716 717 Input Parameters: 718 + linesearch - linesearch context. 719 - snes - The snes instance 720 721 Level: intermediate 722 723 724 .seealso: PetscLineSearchGetSNES(), PetscLineSearchSetVecs() 725 @*/ 726 PetscErrorCode PetscLineSearchSetSNES(PetscLineSearch linesearch, SNES snes){ 727 PetscFunctionBegin; 728 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 729 PetscValidHeaderSpecific(snes,SNES_CLASSID,2); 730 linesearch->snes = snes; 731 PetscFunctionReturn(0); 732 } 733 734 #undef __FUNCT__ 735 #define __FUNCT__ "PetscLineSearchGetSNES" 736 /*@ 737 PetscLineSearchGetSNES - Gets the SNES for the linesearch for function evaluation 738 739 Input Parameters: 740 . linesearch - linesearch context. 741 742 Output Parameters: 743 . snes - The snes instance 744 745 Level: intermediate 746 747 .seealso: PetscLineSearchGetSNES(), PetscLineSearchSetVecs() 748 @*/ 749 PetscErrorCode PetscLineSearchGetSNES(PetscLineSearch linesearch, SNES *snes){ 750 PetscFunctionBegin; 751 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 752 PetscValidPointer(snes, 2); 753 *snes = linesearch->snes; 754 PetscFunctionReturn(0); 755 } 756 757 #undef __FUNCT__ 758 #define __FUNCT__ "PetscLineSearchGetLambda" 759 /*@ 760 PetscLineSearchGetLambda - Gets the last linesearch steplength discovered. 761 762 Input Parameters: 763 . linesearch - linesearch context. 764 765 Output Parameters: 766 . lambda - The last steplength. 767 768 Level: intermediate 769 770 .seealso: PetscLineSearchGetSNES(), PetscLineSearchSetVecs() 771 @*/ 772 PetscErrorCode PetscLineSearchGetLambda(PetscLineSearch linesearch,PetscReal *lambda) 773 { 774 PetscFunctionBegin; 775 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 776 PetscValidPointer(lambda, 2); 777 *lambda = linesearch->lambda; 778 PetscFunctionReturn(0); 779 } 780 781 #undef __FUNCT__ 782 #define __FUNCT__ "PetscLineSearchSetLambda" 783 /*@ 784 PetscLineSearchSetLambda - Sets the linesearch steplength. 785 786 Input Parameters: 787 + linesearch - linesearch context. 788 - lambda - The last steplength. 789 790 Level: intermediate 791 792 .seealso: PetscLineSearchGetLambda() 793 @*/ 794 PetscErrorCode PetscLineSearchSetLambda(PetscLineSearch linesearch, PetscReal lambda) 795 { 796 PetscFunctionBegin; 797 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 798 linesearch->lambda = lambda; 799 PetscFunctionReturn(0); 800 } 801 802 #undef __FUNCT__ 803 #define __FUNCT__ "PetscLineSearchGetTolerances" 804 /*@ 805 PetscLineSearchGetTolerances - Gets the tolerances for the method 806 807 Input Parameters: 808 . linesearch - linesearch context. 809 810 Output Parameters: 811 + steptol - The minimum steplength 812 . rtol - The relative tolerance for iterative line searches 813 . atol - The absolute tolerance for iterative line searches 814 . ltol - The change in lambda tolerance for iterative line searches 815 - max_it - The maximum number of iterations of the line search 816 817 818 Level: advanced 819 820 .seealso: PetscLineSearchSetTolerances() 821 @*/ 822 PetscErrorCode PetscLineSearchGetTolerances(PetscLineSearch linesearch,PetscReal *steptol,PetscReal *maxstep, PetscReal *rtol, PetscReal *atol, PetscReal *ltol, PetscInt *max_its) 823 { 824 PetscFunctionBegin; 825 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 826 if (steptol) { 827 PetscValidPointer(steptol, 2); 828 *steptol = linesearch->steptol; 829 } 830 if (maxstep) { 831 PetscValidPointer(maxstep, 3); 832 *maxstep = linesearch->maxstep; 833 } 834 if (rtol) { 835 PetscValidPointer(rtol, 4); 836 *rtol = linesearch->rtol; 837 } 838 if (atol) { 839 PetscValidPointer(atol, 5); 840 *atol = linesearch->atol; 841 } 842 if (ltol) { 843 PetscValidPointer(ltol, 6); 844 *ltol = linesearch->ltol; 845 } 846 if (max_its) { 847 PetscValidPointer(max_its, 7); 848 *max_its = linesearch->max_its; 849 } 850 PetscFunctionReturn(0); 851 } 852 853 #undef __FUNCT__ 854 #define __FUNCT__ "PetscLineSearchSetTolerances" 855 /*@ 856 PetscLineSearchSetTolerances - Sets the tolerances for the method 857 858 Input Parameters: 859 + linesearch - linesearch context. 860 . steptol - The minimum steplength 861 . rtol - The relative tolerance for iterative line searches 862 . atol - The absolute tolerance for iterative line searches 863 . ltol - The change in lambda tolerance for iterative line searches 864 - max_it - The maximum number of iterations of the line search 865 866 867 Level: advanced 868 869 .seealso: PetscLineSearchGetTolerances() 870 @*/ 871 PetscErrorCode PetscLineSearchSetTolerances(PetscLineSearch linesearch,PetscReal steptol,PetscReal maxstep, PetscReal rtol, PetscReal atol, PetscReal ltol, PetscInt max_its) 872 { 873 PetscFunctionBegin; 874 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 875 linesearch->steptol = steptol; 876 linesearch->maxstep = maxstep; 877 linesearch->rtol = rtol; 878 linesearch->atol = atol; 879 linesearch->ltol = ltol; 880 linesearch->max_its = max_its; 881 PetscFunctionReturn(0); 882 } 883 884 885 #undef __FUNCT__ 886 #define __FUNCT__ "PetscLineSearchGetDamping" 887 /*@ 888 PetscLineSearchGetDamping - Gets the line search damping parameter. 889 890 Input Parameters: 891 . linesearch - linesearch context. 892 893 Output Parameters: 894 . damping - The damping parameter. 895 896 Level: intermediate 897 898 .seealso: PetscLineSearchGetStepTolerance() 899 @*/ 900 901 PetscErrorCode PetscLineSearchGetDamping(PetscLineSearch linesearch,PetscReal *damping) 902 { 903 PetscFunctionBegin; 904 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 905 PetscValidPointer(damping, 2); 906 *damping = linesearch->damping; 907 PetscFunctionReturn(0); 908 } 909 910 #undef __FUNCT__ 911 #define __FUNCT__ "PetscLineSearchSetDamping" 912 /*@ 913 PetscLineSearchSetDamping - Sets the line search damping paramter. 914 915 Input Parameters: 916 . linesearch - linesearch context. 917 . damping - The damping parameter. 918 919 Level: intermediate 920 921 .seealso: PetscLineSearchGetDamping() 922 @*/ 923 PetscErrorCode PetscLineSearchSetDamping(PetscLineSearch linesearch,PetscReal damping) 924 { 925 PetscFunctionBegin; 926 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 927 linesearch->damping = damping; 928 PetscFunctionReturn(0); 929 } 930 931 #undef __FUNCT__ 932 #define __FUNCT__ "PetscLineSearchGetNorms" 933 /*@ 934 PetscLineSearchGetNorms - Gets the norms for for X, Y, and F. 935 936 Input Parameters: 937 . linesearch - linesearch context. 938 939 Output Parameters: 940 + xnorm - The norm of the current solution 941 . fnorm - The norm of the current function 942 - ynorm - The norm of the current update 943 944 Level: intermediate 945 946 .seealso: PetscLineSearchSetNorms() PetscLineSearchGetVecs() 947 @*/ 948 PetscErrorCode PetscLineSearchGetNorms(PetscLineSearch linesearch, PetscReal * xnorm, PetscReal * fnorm, PetscReal * ynorm) 949 { 950 PetscFunctionBegin; 951 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 952 if (xnorm) { 953 *xnorm = linesearch->xnorm; 954 } 955 if (fnorm) { 956 *fnorm = linesearch->fnorm; 957 } 958 if (ynorm) { 959 *ynorm = linesearch->ynorm; 960 } 961 PetscFunctionReturn(0); 962 } 963 964 #undef __FUNCT__ 965 #define __FUNCT__ "PetscLineSearchSetNorms" 966 /*@ 967 PetscLineSearchSetNorms - Gets the computed norms for for X, Y, and F. 968 969 Input Parameters: 970 + linesearch - linesearch context. 971 . xnorm - The norm of the current solution 972 . fnorm - The norm of the current function 973 - ynorm - The norm of the current update 974 975 Level: intermediate 976 977 .seealso: PetscLineSearchGetNorms(), PetscLineSearchSetVecs() 978 @*/ 979 PetscErrorCode PetscLineSearchSetNorms(PetscLineSearch linesearch, PetscReal xnorm, PetscReal fnorm, PetscReal ynorm) 980 { 981 PetscFunctionBegin; 982 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 983 linesearch->xnorm = xnorm; 984 linesearch->fnorm = fnorm; 985 linesearch->ynorm = ynorm; 986 PetscFunctionReturn(0); 987 } 988 989 #undef __FUNCT__ 990 #define __FUNCT__ "PetscLineSearchComputeNorms" 991 /*@ 992 PetscLineSearchComputeNorms - Computes the norms of X, F, and Y. 993 994 Input Parameters: 995 . linesearch - linesearch context. 996 997 Options Database Keys: 998 . -linesearch_norms - turn norm computation on or off. 999 1000 Level: intermediate 1001 1002 .seealso: PetscLineSearchGetNorms, PetscLineSearchSetNorms() 1003 @*/ 1004 PetscErrorCode PetscLineSearchComputeNorms(PetscLineSearch linesearch) 1005 { 1006 PetscErrorCode ierr; 1007 SNES snes; 1008 PetscFunctionBegin; 1009 if (linesearch->norms) { 1010 if (linesearch->ops->vinorm) { 1011 ierr = PetscLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); 1012 ierr = VecNorm(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1013 ierr = VecNorm(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1014 ierr = (*linesearch->ops->vinorm)(snes, linesearch->vec_func, linesearch->vec_sol, &linesearch->fnorm);CHKERRQ(ierr); 1015 } else { 1016 ierr = VecNormBegin(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 1017 ierr = VecNormBegin(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1018 ierr = VecNormBegin(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1019 ierr = VecNormEnd(linesearch->vec_func, NORM_2, &linesearch->fnorm);CHKERRQ(ierr); 1020 ierr = VecNormEnd(linesearch->vec_sol, NORM_2, &linesearch->xnorm);CHKERRQ(ierr); 1021 ierr = VecNormEnd(linesearch->vec_update, NORM_2, &linesearch->ynorm);CHKERRQ(ierr); 1022 } 1023 } 1024 PetscFunctionReturn(0); 1025 } 1026 1027 #undef __FUNCT__ 1028 #define __FUNCT__ "PetscLineSearchGetVecs" 1029 /*@ 1030 PetscLineSearchGetVecs - Gets the vectors from the PetscLineSearch context 1031 1032 Input Parameters: 1033 . linesearch - linesearch context. 1034 1035 Output Parameters: 1036 + X - The old solution 1037 . F - The old function 1038 . Y - The search direction 1039 . W - The new solution 1040 - G - The new function 1041 1042 Level: intermediate 1043 1044 .seealso: PetscLineSearchGetNorms(), PetscLineSearchSetVecs() 1045 @*/ 1046 PetscErrorCode PetscLineSearchGetVecs(PetscLineSearch linesearch,Vec *X,Vec *F, Vec *Y,Vec *W,Vec *G) { 1047 PetscFunctionBegin; 1048 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 1049 if (X) { 1050 PetscValidPointer(X, 2); 1051 *X = linesearch->vec_sol; 1052 } 1053 if (F) { 1054 PetscValidPointer(F, 3); 1055 *F = linesearch->vec_func; 1056 } 1057 if (Y) { 1058 PetscValidPointer(Y, 4); 1059 *Y = linesearch->vec_update; 1060 } 1061 if (W) { 1062 PetscValidPointer(W, 5); 1063 *W = linesearch->vec_sol_new; 1064 } 1065 if (G) { 1066 PetscValidPointer(G, 6); 1067 *G = linesearch->vec_func_new; 1068 } 1069 1070 PetscFunctionReturn(0); 1071 } 1072 1073 #undef __FUNCT__ 1074 #define __FUNCT__ "PetscLineSearchSetVecs" 1075 /*@ 1076 PetscLineSearchSetVecs - Sets the vectors on the PetscLineSearch context 1077 1078 Input Parameters: 1079 + linesearch - linesearch context. 1080 . X - The old solution 1081 . F - The old function 1082 . Y - The search direction 1083 . W - The new solution 1084 - G - The new function 1085 1086 Level: intermediate 1087 1088 .seealso: PetscLineSearchSetNorms(), PetscLineSearchGetVecs() 1089 @*/ 1090 PetscErrorCode PetscLineSearchSetVecs(PetscLineSearch linesearch,Vec X,Vec F,Vec Y,Vec W, Vec G) { 1091 PetscFunctionBegin; 1092 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 1093 if (X) { 1094 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1095 linesearch->vec_sol = X; 1096 } 1097 if (F) { 1098 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 1099 linesearch->vec_func = F; 1100 } 1101 if (Y) { 1102 PetscValidHeaderSpecific(Y,VEC_CLASSID,4); 1103 linesearch->vec_update = Y; 1104 } 1105 if (W) { 1106 PetscValidHeaderSpecific(W,VEC_CLASSID,5); 1107 linesearch->vec_sol_new = W; 1108 } 1109 if (G) { 1110 PetscValidHeaderSpecific(G,VEC_CLASSID,6); 1111 linesearch->vec_func_new = G; 1112 } 1113 1114 PetscFunctionReturn(0); 1115 } 1116 1117 #undef __FUNCT__ 1118 #define __FUNCT__ "PetscLineSearchAppendOptionsPrefix" 1119 /*@C 1120 PetscLineSearchAppendOptionsPrefix - Appends to the prefix used for searching for all 1121 SNES options in the database. 1122 1123 Logically Collective on SNES 1124 1125 Input Parameters: 1126 + snes - the SNES context 1127 - prefix - the prefix to prepend to all option names 1128 1129 Notes: 1130 A hyphen (-) must NOT be given at the beginning of the prefix name. 1131 The first character of all runtime options is AUTOMATICALLY the hyphen. 1132 1133 Level: advanced 1134 1135 .keywords: PetscLineSearch, append, options, prefix, database 1136 1137 .seealso: SNESGetOptionsPrefix() 1138 @*/ 1139 PetscErrorCode PetscLineSearchAppendOptionsPrefix(PetscLineSearch linesearch,const char prefix[]) 1140 { 1141 PetscErrorCode ierr; 1142 1143 PetscFunctionBegin; 1144 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 1145 ierr = PetscObjectAppendOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1146 PetscFunctionReturn(0); 1147 } 1148 1149 #undef __FUNCT__ 1150 #define __FUNCT__ "PetscLineSearchGetOptionsPrefix" 1151 /*@C 1152 PetscLineSearchGetOptionsPrefix - Sets the prefix used for searching for all 1153 PetscLineSearch options in the database. 1154 1155 Not Collective 1156 1157 Input Parameter: 1158 . snes - the SNES context 1159 1160 Output Parameter: 1161 . prefix - pointer to the prefix string used 1162 1163 Notes: On the fortran side, the user should pass in a string 'prefix' of 1164 sufficient length to hold the prefix. 1165 1166 Level: advanced 1167 1168 .keywords: PetscLineSearch, get, options, prefix, database 1169 1170 .seealso: SNESAppendOptionsPrefix() 1171 @*/ 1172 PetscErrorCode PetscLineSearchGetOptionsPrefix(PetscLineSearch linesearch,const char *prefix[]) 1173 { 1174 PetscErrorCode ierr; 1175 1176 PetscFunctionBegin; 1177 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 1178 ierr = PetscObjectGetOptionsPrefix((PetscObject)linesearch,prefix);CHKERRQ(ierr); 1179 PetscFunctionReturn(0); 1180 } 1181 1182 #undef __FUNCT__ 1183 #define __FUNCT__ "PetscLineSearchGetWork" 1184 /*@ 1185 PetscLineSearchGetWork - Gets work vectors for the line search. 1186 1187 Input Parameter: 1188 + linesearch - the PetscLineSearch context 1189 - nwork - the number of work vectors 1190 1191 Level: developer 1192 1193 .keywords: PetscLineSearch, work, vector 1194 1195 .seealso: SNESDefaultGetWork() 1196 @*/ 1197 PetscErrorCode PetscLineSearchGetWork(PetscLineSearch linesearch, PetscInt nwork) 1198 { 1199 PetscErrorCode ierr; 1200 PetscFunctionBegin; 1201 if (linesearch->vec_sol) { 1202 ierr = VecDuplicateVecs(linesearch->vec_sol, nwork, &linesearch->work);CHKERRQ(ierr); 1203 } else { 1204 SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "Cannot get linesearch work-vectors without setting a solution vec!"); 1205 } 1206 PetscFunctionReturn(0); 1207 } 1208 1209 1210 #undef __FUNCT__ 1211 #define __FUNCT__ "PetscLineSearchGetSuccess" 1212 /*@ 1213 PetscLineSearchGetSuccess - Gets the success/failure status of the last line search application 1214 1215 Input Parameters: 1216 . linesearch - linesearch context. 1217 1218 Output Parameters: 1219 . success - The success or failure status. 1220 1221 Level: intermediate 1222 1223 .seealso: PetscLineSearchSetSuccess() 1224 @*/ 1225 PetscErrorCode PetscLineSearchGetSuccess(PetscLineSearch linesearch, PetscBool *success) 1226 { 1227 PetscFunctionBegin; 1228 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 1229 PetscValidPointer(success, 2); 1230 if (success) { 1231 *success = linesearch->success; 1232 } 1233 PetscFunctionReturn(0); 1234 } 1235 1236 #undef __FUNCT__ 1237 #define __FUNCT__ "PetscLineSearchSetSuccess" 1238 /*@ 1239 PetscLineSearchSetSuccess - Sets the success/failure status of the last line search application 1240 1241 Input Parameters: 1242 + linesearch - linesearch context. 1243 - success - The success or failure status. 1244 1245 Level: intermediate 1246 1247 .seealso: PetscLineSearchGetSuccess() 1248 @*/ 1249 PetscErrorCode PetscLineSearchSetSuccess(PetscLineSearch linesearch, PetscBool success) 1250 { 1251 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 1252 PetscFunctionBegin; 1253 linesearch->success = success; 1254 PetscFunctionReturn(0); 1255 } 1256 1257 #undef __FUNCT__ 1258 #define __FUNCT__ "PetscLineSearchSetVIFunctions" 1259 /*@C 1260 PetscLineSearchSetVIFunctions - Sets VI-specific functions for line search computation. 1261 1262 Input Parameters: 1263 + snes - nonlinear context obtained from SNESCreate() 1264 . projectfunc - function for projecting the function to the bounds 1265 - normfunc - function for computing the norm of an active set 1266 1267 Logically Collective on SNES 1268 1269 Calling sequence of projectfunc: 1270 .vb 1271 projectfunc (SNES snes, Vec X) 1272 .ve 1273 1274 Input parameters for projectfunc: 1275 + snes - nonlinear context 1276 - X - current solution 1277 1278 Output parameters for func: 1279 . X - Projected solution 1280 1281 Calling sequence of normfunc: 1282 .vb 1283 projectfunc (SNES snes, Vec X, Vec F, PetscScalar * fnorm) 1284 .ve 1285 1286 Input parameters for projectfunc: 1287 + snes - nonlinear context 1288 . X - current solution 1289 - F - current residual 1290 1291 Output parameters for func: 1292 . fnorm - VI-specific norm of the function 1293 1294 1295 Level: developer 1296 1297 .keywords: SNES, line search, VI, nonlinear, set, line search 1298 1299 .seealso: PetscLineSearchGetVIFunctions(), PetscLineSearchSetPostCheck(), PetscLineSearchSetPreCheck() 1300 @*/ 1301 extern PetscErrorCode PetscLineSearchSetVIFunctions(PetscLineSearch linesearch, PetscLineSearchVIProjectFunc projectfunc, PetscLineSearchVINormFunc normfunc) 1302 { 1303 PetscFunctionBegin; 1304 PetscValidHeaderSpecific(linesearch,PETSCLINESEARCH_CLASSID,1); 1305 if (projectfunc) linesearch->ops->viproject = projectfunc; 1306 if (normfunc) linesearch->ops->vinorm = normfunc; 1307 PetscFunctionReturn(0); 1308 } 1309 1310 #undef __FUNCT__ 1311 #define __FUNCT__ "PetscLineSearchGetVIFunctions" 1312 /*@C 1313 PetscLineSearchGetVIFunctions - Sets VI-specific functions for line search computation. 1314 1315 Input Parameters: 1316 . snes - nonlinear context obtained from SNESCreate() 1317 1318 Output Parameters: 1319 + projectfunc - function for projecting the function to the bounds 1320 - normfunc - function for computing the norm of an active set 1321 1322 Logically Collective on SNES 1323 1324 Level: developer 1325 1326 .keywords: SNES, line search, VI, nonlinear, get, line search 1327 1328 .seealso: PetscLineSearchSetVIFunctions(), PetscLineSearchGetPostCheck(), PetscLineSearchGetPreCheck() 1329 @*/ 1330 extern PetscErrorCode PetscLineSearchGetVIFunctions(PetscLineSearch linesearch, PetscLineSearchVIProjectFunc *projectfunc, PetscLineSearchVINormFunc *normfunc) 1331 { 1332 PetscFunctionBegin; 1333 if (projectfunc) *projectfunc = linesearch->ops->viproject; 1334 if (normfunc) *normfunc = linesearch->ops->vinorm; 1335 PetscFunctionReturn(0); 1336 } 1337 1338 #undef __FUNCT__ 1339 #define __FUNCT__ "PetscLineSearchRegister" 1340 /*@C 1341 PetscLineSearchRegister - See PetscLineSearchRegisterDynamic() 1342 1343 Level: advanced 1344 @*/ 1345 PetscErrorCode PetscLineSearchRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(PetscLineSearch)) 1346 { 1347 char fullname[PETSC_MAX_PATH_LEN]; 1348 PetscErrorCode ierr; 1349 1350 PetscFunctionBegin; 1351 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 1352 ierr = PetscFListAdd(&PetscLineSearchList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 1353 PetscFunctionReturn(0); 1354 } 1355