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