1 /*$Id: ls.c,v 1.156 2000/07/23 16:42:19 bsmith Exp bsmith $*/ 2 3 #include "src/snes/impls/ls/ls.h" 4 5 /* 6 Checks if J^T F = 0 which implies we've found a local minimum of the function, 7 but not a zero. In the case when one cannot compute J^T F we use the fact that 8 0 = (J^T F)^T W = F^T J W iff W not in the null space of J 9 */ 10 #undef __FUNC__ 11 #define __FUNC__ /*<a name="SNESLSCheckLocalMin_Private"></a>*/"SNESLSCheckLocalMin_Private" 12 int SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin) 13 { 14 PetscReal a1; 15 int ierr; 16 17 PetscFunctionBegin; 18 *ismin = PETSC_FALSE; 19 /* Compute || J^T F|| */ 20 ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 21 ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 22 PLogInfo(0,"SNESSolve_EQ_LS: || J^T F|| %g || F || %g ||J^T F||/|| F || %g\n",a1,fnorm,a1/fnorm); 23 if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 24 PetscFunctionReturn(0); 25 } 26 27 /* -------------------------------------------------------------------- 28 29 This file implements a truncated Newton method with a line search, 30 for solving a system of nonlinear equations, using the SLES, Vec, 31 and Mat interfaces for linear solvers, vectors, and matrices, 32 respectively. 33 34 The following basic routines are required for each nonlinear solver: 35 SNESCreate_XXX() - Creates a nonlinear solver context 36 SNESSetFromOptions_XXX() - Sets runtime options 37 SNESSolve_XXX() - Solves the nonlinear system 38 SNESDestroy_XXX() - Destroys the nonlinear solver context 39 The suffix "_XXX" denotes a particular implementation, in this case 40 we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving 41 systems of nonlinear equations (EQ) with a line search (LS) method. 42 These routines are actually called via the common user interface 43 routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 44 SNESDestroy(), so the application code interface remains identical 45 for all nonlinear solvers. 46 47 Another key routine is: 48 SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 49 by setting data structures and options. The interface routine SNESSetUp() 50 is not usually called directly by the user, but instead is called by 51 SNESSolve() if necessary. 52 53 Additional basic routines are: 54 SNESPrintHelp_XXX() - Prints nonlinear solver runtime options 55 SNESView_XXX() - Prints details of runtime options that 56 have actually been used. 57 These are called by application codes via the interface routines 58 SNESPrintHelp() and SNESView(). 59 60 The various types of solvers (preconditioners, Krylov subspace methods, 61 nonlinear solvers, timesteppers) are all organized similarly, so the 62 above description applies to these categories also. 63 64 -------------------------------------------------------------------- */ 65 /* 66 SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton 67 method with a line search. 68 69 Input Parameters: 70 . snes - the SNES context 71 72 Output Parameter: 73 . outits - number of iterations until termination 74 75 Application Interface Routine: SNESSolve() 76 77 Notes: 78 This implements essentially a truncated Newton method with a 79 line search. By default a cubic backtracking line search 80 is employed, as described in the text "Numerical Methods for 81 Unconstrained Optimization and Nonlinear Equations" by Dennis 82 and Schnabel. 83 */ 84 #undef __FUNC__ 85 #define __FUNC__ /*<a name=""></a>*/"SNESSolve_EQ_LS" 86 int SNESSolve_EQ_LS(SNES snes,int *outits) 87 { 88 SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 89 int maxits,i,ierr,lits,lsfail; 90 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 91 PetscReal fnorm,gnorm,xnorm,ynorm; 92 Vec Y,X,F,G,W,TMP; 93 94 PetscFunctionBegin; 95 snes->reason = SNES_CONVERGED_ITERATING; 96 97 maxits = snes->max_its; /* maximum number of iterations */ 98 X = snes->vec_sol; /* solution vector */ 99 F = snes->vec_func; /* residual vector */ 100 Y = snes->work[0]; /* work vectors */ 101 G = snes->work[1]; 102 W = snes->work[2]; 103 104 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 105 snes->iter = 0; 106 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 107 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 108 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 109 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 110 snes->norm = fnorm; 111 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 112 SNESLogConvHistory(snes,fnorm,0); 113 SNESMonitor(snes,0,fnorm); 114 115 if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);} 116 117 /* set parameter for default relative tolerance convergence test */ 118 snes->ttol = fnorm*snes->rtol; 119 120 for (i=0; i<maxits; i++) { 121 122 /* Solve J Y = F, where J is Jacobian matrix */ 123 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 124 ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 125 ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr); 126 /* should check what happened to the linear solve? */ 127 snes->linear_its += lits; 128 PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits); 129 130 /* Compute a (scaled) negative update in the line search routine: 131 Y <- X - lambda*Y 132 and evaluate G(Y) = function(Y)) 133 */ 134 ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 135 ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr); 136 PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail); 137 if (lsfail) { 138 PetscTruth ismin; 139 snes->nfailures++; 140 snes->reason = SNES_DIVERGED_LS_FAILURE; 141 ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); 142 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 143 break; 144 } 145 146 TMP = F; F = G; snes->vec_func_always = F; G = TMP; 147 TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 148 fnorm = gnorm; 149 150 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 151 snes->iter = i+1; 152 snes->norm = fnorm; 153 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 154 SNESLogConvHistory(snes,fnorm,lits); 155 SNESMonitor(snes,i+1,fnorm); 156 157 /* Test for convergence */ 158 if (snes->converged) { 159 ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 160 ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 161 if (snes->reason) { 162 break; 163 } 164 } 165 } 166 if (X != snes->vec_sol) { 167 ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 168 } 169 if (F != snes->vec_func) { 170 ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 171 } 172 snes->vec_sol_always = snes->vec_sol; 173 snes->vec_func_always = snes->vec_func; 174 if (i == maxits) { 175 PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits); 176 i--; 177 snes->reason = SNES_DIVERGED_MAX_IT; 178 } 179 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 180 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 181 *outits = i+1; 182 PetscFunctionReturn(0); 183 } 184 /* -------------------------------------------------------------------------- */ 185 /* 186 SNESSetUp_EQ_LS - Sets up the internal data structures for the later use 187 of the SNESEQLS nonlinear solver. 188 189 Input Parameter: 190 . snes - the SNES context 191 . x - the solution vector 192 193 Application Interface Routine: SNESSetUp() 194 195 Notes: 196 For basic use of the SNES solvers the user need not explicitly call 197 SNESSetUp(), since these actions will automatically occur during 198 the call to SNESSolve(). 199 */ 200 #undef __FUNC__ 201 #define __FUNC__ /*<a name=""></a>*/"SNESSetUp_EQ_LS" 202 int SNESSetUp_EQ_LS(SNES snes) 203 { 204 int ierr; 205 206 PetscFunctionBegin; 207 snes->nwork = 4; 208 ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 209 PLogObjectParents(snes,snes->nwork,snes->work); 210 snes->vec_sol_update_always = snes->work[3]; 211 PetscFunctionReturn(0); 212 } 213 /* -------------------------------------------------------------------------- */ 214 /* 215 SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created 216 with SNESCreate_EQ_LS(). 217 218 Input Parameter: 219 . snes - the SNES context 220 221 Application Interface Routine: SNESDestroy() 222 */ 223 #undef __FUNC__ 224 #define __FUNC__ /*<a name=""></a>*/"SNESDestroy_EQ_LS" 225 int SNESDestroy_EQ_LS(SNES snes) 226 { 227 int ierr; 228 229 PetscFunctionBegin; 230 if (snes->nwork) { 231 ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 232 } 233 ierr = PetscFree(snes->data);CHKERRQ(ierr); 234 PetscFunctionReturn(0); 235 } 236 /* -------------------------------------------------------------------------- */ 237 #undef __FUNC__ 238 #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearch" 239 240 /*@C 241 SNESNoLineSearch - This routine is not a line search at all; 242 it simply uses the full Newton step. Thus, this routine is intended 243 to serve as a template and is not recommended for general use. 244 245 Collective on SNES and Vec 246 247 Input Parameters: 248 + snes - nonlinear context 249 . lsctx - optional context for line search (not used here) 250 . x - current iterate 251 . f - residual evaluated at x 252 . y - search direction (contains new iterate on output) 253 . w - work vector 254 - fnorm - 2-norm of f 255 256 Output Parameters: 257 + g - residual evaluated at new iterate y 258 . y - new iterate (contains search direction on input) 259 . gnorm - 2-norm of g 260 . ynorm - 2-norm of search length 261 - flag - set to 0, indicating a successful line search 262 263 Options Database Key: 264 . -snes_eq_ls basic - Activates SNESNoLineSearch() 265 266 Level: advanced 267 268 .keywords: SNES, nonlinear, line search, cubic 269 270 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 271 SNESSetLineSearch(), SNESNoLineSearchNoNorms() 272 @*/ 273 int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 274 PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 275 { 276 int ierr; 277 Scalar mone = -1.0; 278 SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 279 PetscTruth change_y = PETSC_FALSE; 280 281 PetscFunctionBegin; 282 *flag = 0; 283 PLogEventBegin(SNES_LineSearch,snes,x,f,g); 284 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 285 ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 286 if (neP->CheckStep) { 287 ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 288 } 289 ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 290 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 291 PLogEventEnd(SNES_LineSearch,snes,x,f,g); 292 PetscFunctionReturn(0); 293 } 294 /* -------------------------------------------------------------------------- */ 295 #undef __FUNC__ 296 #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearchNoNorms" 297 298 /*@C 299 SNESNoLineSearchNoNorms - This routine is not a line search at 300 all; it simply uses the full Newton step. This version does not 301 even compute the norm of the function or search direction; this 302 is intended only when you know the full step is fine and are 303 not checking for convergence of the nonlinear iteration (for 304 example, you are running always for a fixed number of Newton steps). 305 306 Collective on SNES and Vec 307 308 Input Parameters: 309 + snes - nonlinear context 310 . lsctx - optional context for line search (not used here) 311 . x - current iterate 312 . f - residual evaluated at x 313 . y - search direction (contains new iterate on output) 314 . w - work vector 315 - fnorm - 2-norm of f 316 317 Output Parameters: 318 + g - residual evaluated at new iterate y 319 . gnorm - not changed 320 . ynorm - not changed 321 - flag - set to 0, indicating a successful line search 322 323 Options Database Key: 324 . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 325 326 Notes: 327 SNESNoLineSearchNoNorms() must be used in conjunction with 328 either the options 329 $ -snes_no_convergence_test -snes_max_it <its> 330 or alternatively a user-defined custom test set via 331 SNESSetConvergenceTest(); or a -snes_max_it of 1, 332 otherwise, the SNES solver will generate an error. 333 334 During the final iteration this will not evaluate the function at 335 the solution point. This is to save a function evaluation while 336 using pseudo-timestepping. 337 338 The residual norms printed by monitoring routines such as 339 SNESDefaultMonitor() (as activated via -snes_monitor) will not be 340 correct, since they are not computed. 341 342 Level: advanced 343 344 .keywords: SNES, nonlinear, line search, cubic 345 346 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 347 SNESSetLineSearch(), SNESNoLineSearch() 348 @*/ 349 int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 350 PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 351 { 352 int ierr; 353 Scalar mone = -1.0; 354 SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 355 PetscTruth change_y = PETSC_FALSE; 356 357 PetscFunctionBegin; 358 *flag = 0; 359 PLogEventBegin(SNES_LineSearch,snes,x,f,g); 360 ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 361 if (neP->CheckStep) { 362 ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 363 } 364 365 /* don't evaluate function the last time through */ 366 if (snes->iter < snes->max_its-1) { 367 ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 368 } 369 PLogEventEnd(SNES_LineSearch,snes,x,f,g); 370 PetscFunctionReturn(0); 371 } 372 /* -------------------------------------------------------------------------- */ 373 #undef __FUNC__ 374 #define __FUNC__ /*<a name=""></a>*/"SNESCubicLineSearch" 375 /*@C 376 SNESCubicLineSearch - Performs a cubic line search (default line search method). 377 378 Collective on SNES 379 380 Input Parameters: 381 + snes - nonlinear context 382 . lsctx - optional context for line search (not used here) 383 . x - current iterate 384 . f - residual evaluated at x 385 . y - search direction (contains new iterate on output) 386 . w - work vector 387 - fnorm - 2-norm of f 388 389 Output Parameters: 390 + g - residual evaluated at new iterate y 391 . y - new iterate (contains search direction on input) 392 . gnorm - 2-norm of g 393 . ynorm - 2-norm of search length 394 - flag - 0 if line search succeeds; -1 on failure. 395 396 Options Database Key: 397 $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 398 399 Notes: 400 This line search is taken from "Numerical Methods for Unconstrained 401 Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 402 403 Level: advanced 404 405 .keywords: SNES, nonlinear, line search, cubic 406 407 .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 408 @*/ 409 int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 410 PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 411 { 412 /* 413 Note that for line search purposes we work with with the related 414 minimization problem: 415 min z(x): R^n -> R, 416 where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 417 */ 418 419 PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2; 420 PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 421 #if defined(PETSC_USE_COMPLEX) 422 Scalar cinitslope,clambda; 423 #endif 424 int ierr,count; 425 SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 426 Scalar mone = -1.0,scale; 427 PetscTruth change_y = PETSC_FALSE; 428 429 PetscFunctionBegin; 430 PLogEventBegin(SNES_LineSearch,snes,x,f,g); 431 *flag = 0; 432 alpha = neP->alpha; 433 maxstep = neP->maxstep; 434 steptol = neP->steptol; 435 436 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 437 if (*ynorm < snes->atol) { 438 PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 439 *gnorm = fnorm; 440 ierr = VecCopy(x,y);CHKERRQ(ierr); 441 ierr = VecCopy(f,g);CHKERRQ(ierr); 442 goto theend1; 443 } 444 if (*ynorm > maxstep) { /* Step too big, so scale back */ 445 scale = maxstep/(*ynorm); 446 #if defined(PETSC_USE_COMPLEX) 447 PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm); 448 #else 449 PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm); 450 #endif 451 ierr = VecScale(&scale,y);CHKERRQ(ierr); 452 *ynorm = maxstep; 453 } 454 minlambda = steptol/(*ynorm); 455 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 456 #if defined(PETSC_USE_COMPLEX) 457 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 458 initslope = PetscRealPart(cinitslope); 459 #else 460 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 461 #endif 462 if (initslope > 0.0) initslope = -initslope; 463 if (initslope == 0.0) initslope = -1.0; 464 465 ierr = VecCopy(y,w);CHKERRQ(ierr); 466 ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 467 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 468 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 469 if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 470 ierr = VecCopy(w,y);CHKERRQ(ierr); 471 PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 472 goto theend1; 473 } 474 475 /* Fit points with quadratic */ 476 lambda = 1.0; 477 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 478 lambdaprev = lambda; 479 gnormprev = *gnorm; 480 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 481 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 482 else lambda = lambdatemp; 483 ierr = VecCopy(x,w);CHKERRQ(ierr); 484 lambdaneg = -lambda; 485 #if defined(PETSC_USE_COMPLEX) 486 clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 487 #else 488 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 489 #endif 490 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 491 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 492 if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 493 ierr = VecCopy(w,y);CHKERRQ(ierr); 494 PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 495 goto theend1; 496 } 497 498 /* Fit points with cubic */ 499 count = 1; 500 while (1) { 501 if (lambda <= minlambda) { /* bad luck; use full step */ 502 PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 503 PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 504 ierr = VecCopy(w,y);CHKERRQ(ierr); 505 *flag = -1; break; 506 } 507 t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope; 508 t2 = (gnormprev*gnormprev - fnorm*fnorm)*0.5 - lambdaprev*initslope; 509 a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 510 b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 511 d = b*b - 3*a*initslope; 512 if (d < 0.0) d = 0.0; 513 if (a == 0.0) { 514 lambdatemp = -initslope/(2.0*b); 515 } else { 516 lambdatemp = (-b + sqrt(d))/(3.0*a); 517 } 518 lambdaprev = lambda; 519 gnormprev = *gnorm; 520 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 521 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 522 else lambda = lambdatemp; 523 ierr = VecCopy(x,w);CHKERRQ(ierr); 524 lambdaneg = -lambda; 525 #if defined(PETSC_USE_COMPLEX) 526 clambda = lambdaneg; 527 ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 528 #else 529 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 530 #endif 531 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 532 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 533 if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */ 534 ierr = VecCopy(w,y);CHKERRQ(ierr); 535 PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 536 goto theend1; 537 } else { 538 PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 539 } 540 count++; 541 } 542 theend1: 543 /* Optional user-defined check for line search step validity */ 544 if (neP->CheckStep) { 545 ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 546 if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 547 ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 548 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 549 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 550 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 551 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 552 } 553 } 554 PLogEventEnd(SNES_LineSearch,snes,x,f,g); 555 PetscFunctionReturn(0); 556 } 557 /* -------------------------------------------------------------------------- */ 558 #undef __FUNC__ 559 #define __FUNC__ /*<a name=""></a>*/"SNESQuadraticLineSearch" 560 /*@C 561 SNESQuadraticLineSearch - Performs a quadratic line search. 562 563 Collective on SNES and Vec 564 565 Input Parameters: 566 + snes - the SNES context 567 . lsctx - optional context for line search (not used here) 568 . x - current iterate 569 . f - residual evaluated at x 570 . y - search direction (contains new iterate on output) 571 . w - work vector 572 - fnorm - 2-norm of f 573 574 Output Parameters: 575 + g - residual evaluated at new iterate y 576 . y - new iterate (contains search direction on input) 577 . gnorm - 2-norm of g 578 . ynorm - 2-norm of search length 579 - flag - 0 if line search succeeds; -1 on failure. 580 581 Options Database Key: 582 . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 583 584 Notes: 585 Use SNESSetLineSearch() to set this routine within the SNESEQLS method. 586 587 Level: advanced 588 589 .keywords: SNES, nonlinear, quadratic, line search 590 591 .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 592 @*/ 593 int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 594 PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 595 { 596 /* 597 Note that for line search purposes we work with with the related 598 minimization problem: 599 min z(x): R^n -> R, 600 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 601 */ 602 PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 603 #if defined(PETSC_USE_COMPLEX) 604 Scalar cinitslope,clambda; 605 #endif 606 int ierr,count; 607 SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 608 Scalar mone = -1.0,scale; 609 PetscTruth change_y = PETSC_FALSE; 610 611 PetscFunctionBegin; 612 PLogEventBegin(SNES_LineSearch,snes,x,f,g); 613 *flag = 0; 614 alpha = neP->alpha; 615 maxstep = neP->maxstep; 616 steptol = neP->steptol; 617 618 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 619 if (*ynorm < snes->atol) { 620 PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 621 *gnorm = fnorm; 622 ierr = VecCopy(x,y);CHKERRQ(ierr); 623 ierr = VecCopy(f,g);CHKERRQ(ierr); 624 goto theend2; 625 } 626 if (*ynorm > maxstep) { /* Step too big, so scale back */ 627 scale = maxstep/(*ynorm); 628 ierr = VecScale(&scale,y);CHKERRQ(ierr); 629 *ynorm = maxstep; 630 } 631 minlambda = steptol/(*ynorm); 632 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 633 #if defined(PETSC_USE_COMPLEX) 634 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 635 initslope = PetscRealPart(cinitslope); 636 #else 637 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 638 #endif 639 if (initslope > 0.0) initslope = -initslope; 640 if (initslope == 0.0) initslope = -1.0; 641 642 ierr = VecCopy(y,w);CHKERRQ(ierr); 643 ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 644 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 645 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 646 if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 647 ierr = VecCopy(w,y);CHKERRQ(ierr); 648 PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 649 goto theend2; 650 } 651 652 /* Fit points with quadratic */ 653 lambda = 1.0; 654 count = 1; 655 while (1) { 656 if (lambda <= minlambda) { /* bad luck; use full step */ 657 PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 658 PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 659 fnorm,*gnorm,*ynorm,lambda,initslope); 660 ierr = VecCopy(w,y);CHKERRQ(ierr); 661 *flag = -1; break; 662 } 663 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 664 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 665 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 666 else lambda = lambdatemp; 667 ierr = VecCopy(x,w);CHKERRQ(ierr); 668 lambdaneg = -lambda; 669 #if defined(PETSC_USE_COMPLEX) 670 clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 671 #else 672 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 673 #endif 674 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 675 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 676 if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 677 ierr = VecCopy(w,y);CHKERRQ(ierr); 678 PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 679 break; 680 } 681 count++; 682 } 683 theend2: 684 /* Optional user-defined check for line search step validity */ 685 if (neP->CheckStep) { 686 ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 687 if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 688 ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 689 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 690 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 691 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 692 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 693 } 694 } 695 PLogEventEnd(SNES_LineSearch,snes,x,f,g); 696 PetscFunctionReturn(0); 697 } 698 /* -------------------------------------------------------------------------- */ 699 #undef __FUNC__ 700 #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch" 701 /*@C 702 SNESSetLineSearch - Sets the line search routine to be used 703 by the method SNESEQLS. 704 705 Input Parameters: 706 + snes - nonlinear context obtained from SNESCreate() 707 . lsctx - optional user-defined context for use by line search 708 - func - pointer to int function 709 710 Collective on SNES 711 712 Available Routines: 713 + SNESCubicLineSearch() - default line search 714 . SNESQuadraticLineSearch() - quadratic line search 715 . SNESNoLineSearch() - the full Newton step (actually not a line search) 716 - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 717 718 Options Database Keys: 719 + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 720 . -snes_eq_ls_alpha <alpha> - Sets alpha 721 . -snes_eq_ls_maxstep <max> - Sets maxstep 722 - -snes_eq_ls_steptol <steptol> - Sets steptol 723 724 Calling sequence of func: 725 .vb 726 func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 727 PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag) 728 .ve 729 730 Input parameters for func: 731 + snes - nonlinear context 732 . lsctx - optional user-defined context for line search 733 . x - current iterate 734 . f - residual evaluated at x 735 . y - search direction (contains new iterate on output) 736 . w - work vector 737 - fnorm - 2-norm of f 738 739 Output parameters for func: 740 + g - residual evaluated at new iterate y 741 . y - new iterate (contains search direction on input) 742 . gnorm - 2-norm of g 743 . ynorm - 2-norm of search length 744 - flag - set to 0 if the line search succeeds; a nonzero integer 745 on failure. 746 747 Level: advanced 748 749 .keywords: SNES, nonlinear, set, line search, routine 750 751 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 752 @*/ 753 int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 754 { 755 int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*); 756 757 PetscFunctionBegin; 758 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 759 if (f) { 760 ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 761 } 762 PetscFunctionReturn(0); 763 } 764 /* -------------------------------------------------------------------------- */ 765 EXTERN_C_BEGIN 766 #undef __FUNC__ 767 #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch_LS" 768 int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 769 double,double*,double*,int*),void *lsctx) 770 { 771 PetscFunctionBegin; 772 ((SNES_EQ_LS *)(snes->data))->LineSearch = func; 773 ((SNES_EQ_LS *)(snes->data))->lsP = lsctx; 774 PetscFunctionReturn(0); 775 } 776 EXTERN_C_END 777 /* -------------------------------------------------------------------------- */ 778 #undef __FUNC__ 779 #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck" 780 /*@C 781 SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 782 by the line search routine in the Newton-based method SNESEQLS. 783 784 Input Parameters: 785 + snes - nonlinear context obtained from SNESCreate() 786 . func - pointer to int function 787 - checkctx - optional user-defined context for use by step checking routine 788 789 Collective on SNES 790 791 Calling sequence of func: 792 .vb 793 int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 794 .ve 795 where func returns an error code of 0 on success and a nonzero 796 on failure. 797 798 Input parameters for func: 799 + snes - nonlinear context 800 . checkctx - optional user-defined context for use by step checking routine 801 - x - current candidate iterate 802 803 Output parameters for func: 804 + x - current iterate (possibly modified) 805 - flag - flag indicating whether x has been modified (either 806 PETSC_TRUE of PETSC_FALSE) 807 808 Level: advanced 809 810 Notes: 811 SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 812 iterate computed by the line search checking routine. In particular, 813 these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 814 to the checking routine, and then (3) compute the corresponding nonlinear 815 function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 816 817 SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 818 new iterate computed by the line search checking routine. In particular, 819 these routines (1) compute a candidate iterate u_{i+1} as well as a 820 candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 821 routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 822 were made to the candidate iterate in the checking routine (as indicated 823 by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 824 very costly, so use this feature with caution! 825 826 .keywords: SNES, nonlinear, set, line search check, step check, routine 827 828 .seealso: SNESSetLineSearch() 829 @*/ 830 int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 831 { 832 int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 833 834 PetscFunctionBegin; 835 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr); 836 if (f) { 837 ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 838 } 839 PetscFunctionReturn(0); 840 } 841 /* -------------------------------------------------------------------------- */ 842 EXTERN_C_BEGIN 843 #undef __FUNC__ 844 #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck_LS" 845 int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 846 { 847 PetscFunctionBegin; 848 ((SNES_EQ_LS *)(snes->data))->CheckStep = func; 849 ((SNES_EQ_LS *)(snes->data))->checkP = checkctx; 850 PetscFunctionReturn(0); 851 } 852 EXTERN_C_END 853 /* -------------------------------------------------------------------------- */ 854 /* 855 SNESPrintHelp_EQ_LS - Prints all options for the SNESEQLS method. 856 857 Input Parameter: 858 . snes - the SNES context 859 860 Application Interface Routine: SNESPrintHelp() 861 */ 862 #undef __FUNC__ 863 #define __FUNC__ /*<a name=""></a>*/"SNESPrintHelp_EQ_LS" 864 static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 865 { 866 SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 867 int ierr; 868 869 PetscFunctionBegin; 870 ierr = (*PetscHelpPrintf)(snes->comm," method SNESEQLS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr); 871 ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr); 872 ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr); 873 ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr); 874 ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr); 875 PetscFunctionReturn(0); 876 } 877 /* -------------------------------------------------------------------------- */ 878 /* 879 SNESView_EQ_LS - Prints info from the SNESEQLS data structure. 880 881 Input Parameters: 882 . SNES - the SNES context 883 . viewer - visualization context 884 885 Application Interface Routine: SNESView() 886 */ 887 #undef __FUNC__ 888 #define __FUNC__ /*<a name=""></a>*/"SNESView_EQ_LS" 889 static int SNESView_EQ_LS(SNES snes,Viewer viewer) 890 { 891 SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 892 char *cstr; 893 int ierr; 894 PetscTruth isascii; 895 896 PetscFunctionBegin; 897 ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 898 if (isascii) { 899 if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 900 else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 901 else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 902 else cstr = "unknown"; 903 ierr = ViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 904 ierr = ViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 905 } else { 906 SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 907 } 908 PetscFunctionReturn(0); 909 } 910 /* -------------------------------------------------------------------------- */ 911 /* 912 SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method. 913 914 Input Parameter: 915 . snes - the SNES context 916 917 Application Interface Routine: SNESSetFromOptions() 918 */ 919 #undef __FUNC__ 920 #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions_EQ_LS" 921 static int SNESSetFromOptions_EQ_LS(SNES snes) 922 { 923 SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 924 char ver[16]; 925 double tmp; 926 int ierr; 927 PetscTruth flg; 928 929 PetscFunctionBegin; 930 ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp,&flg);CHKERRQ(ierr); 931 if (flg) { 932 ls->alpha = tmp; 933 } 934 ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp,&flg);CHKERRQ(ierr); 935 if (flg) { 936 ls->maxstep = tmp; 937 } 938 ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp,&flg);CHKERRQ(ierr); 939 if (flg) { 940 ls->steptol = tmp; 941 } 942 ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16,&flg);CHKERRQ(ierr); 943 if (flg) { 944 PetscTruth isbasic,isnonorms,isquad,iscubic; 945 946 ierr = PetscStrcmp(ver,"basic",&isbasic);CHKERRQ(ierr); 947 ierr = PetscStrcmp(ver,"basicnonorms",&isnonorms);CHKERRQ(ierr); 948 ierr = PetscStrcmp(ver,"quadratic",&isquad);CHKERRQ(ierr); 949 ierr = PetscStrcmp(ver,"cubic",&iscubic);CHKERRQ(ierr); 950 951 if (isbasic) { 952 ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 953 } else if (isnonorms) { 954 ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 955 } else if (isquad) { 956 ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 957 } else if (iscubic) { 958 ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 959 } 960 else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 961 } 962 PetscFunctionReturn(0); 963 } 964 /* -------------------------------------------------------------------------- */ 965 /* 966 SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method, 967 SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver 968 context, SNES, that was created within SNESCreate(). 969 970 971 Input Parameter: 972 . snes - the SNES context 973 974 Application Interface Routine: SNESCreate() 975 */ 976 EXTERN_C_BEGIN 977 #undef __FUNC__ 978 #define __FUNC__ /*<a name=""></a>*/"SNESCreate_EQ_LS" 979 int SNESCreate_EQ_LS(SNES snes) 980 { 981 int ierr; 982 SNES_EQ_LS *neP; 983 984 PetscFunctionBegin; 985 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 986 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 987 } 988 989 snes->setup = SNESSetUp_EQ_LS; 990 snes->solve = SNESSolve_EQ_LS; 991 snes->destroy = SNESDestroy_EQ_LS; 992 snes->converged = SNESConverged_EQ_LS; 993 snes->printhelp = SNESPrintHelp_EQ_LS; 994 snes->setfromoptions = SNESSetFromOptions_EQ_LS; 995 snes->view = SNESView_EQ_LS; 996 snes->nwork = 0; 997 998 neP = PetscNew(SNES_EQ_LS);CHKPTRQ(neP); 999 PLogObjectMemory(snes,sizeof(SNES_EQ_LS)); 1000 snes->data = (void*)neP; 1001 neP->alpha = 1.e-4; 1002 neP->maxstep = 1.e8; 1003 neP->steptol = 1.e-12; 1004 neP->LineSearch = SNESCubicLineSearch; 1005 neP->lsP = PETSC_NULL; 1006 neP->CheckStep = PETSC_NULL; 1007 neP->checkP = PETSC_NULL; 1008 1009 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS", 1010 SNESSetLineSearch_LS);CHKERRQ(ierr); 1011 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS", 1012 SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 1013 1014 PetscFunctionReturn(0); 1015 } 1016 EXTERN_C_END 1017 1018 1019 1020