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