1 /*$Id: ls.c,v 1.158 2000/07/24 03:50:34 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 near zero implies found a local minimum\n",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 near zero implies found a local minimum\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,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 293 { 294 int ierr; 295 Scalar mone = -1.0; 296 SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 297 PetscTruth change_y = PETSC_FALSE; 298 299 PetscFunctionBegin; 300 *flag = 0; 301 PLogEventBegin(SNES_LineSearch,snes,x,f,g); 302 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 303 ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 304 if (neP->CheckStep) { 305 ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 306 } 307 ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 308 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 309 PLogEventEnd(SNES_LineSearch,snes,x,f,g); 310 PetscFunctionReturn(0); 311 } 312 /* -------------------------------------------------------------------------- */ 313 #undef __FUNC__ 314 #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearchNoNorms" 315 316 /*@C 317 SNESNoLineSearchNoNorms - This routine is not a line search at 318 all; it simply uses the full Newton step. This version does not 319 even compute the norm of the function or search direction; this 320 is intended only when you know the full step is fine and are 321 not checking for convergence of the nonlinear iteration (for 322 example, you are running always for a fixed number of Newton steps). 323 324 Collective on SNES and Vec 325 326 Input Parameters: 327 + snes - nonlinear context 328 . lsctx - optional context for line search (not used here) 329 . x - current iterate 330 . f - residual evaluated at x 331 . y - search direction (contains new iterate on output) 332 . w - work vector 333 - fnorm - 2-norm of f 334 335 Output Parameters: 336 + g - residual evaluated at new iterate y 337 . gnorm - not changed 338 . ynorm - not changed 339 - flag - set to 0, indicating a successful line search 340 341 Options Database Key: 342 . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 343 344 Notes: 345 SNESNoLineSearchNoNorms() must be used in conjunction with 346 either the options 347 $ -snes_no_convergence_test -snes_max_it <its> 348 or alternatively a user-defined custom test set via 349 SNESSetConvergenceTest(); or a -snes_max_it of 1, 350 otherwise, the SNES solver will generate an error. 351 352 During the final iteration this will not evaluate the function at 353 the solution point. This is to save a function evaluation while 354 using pseudo-timestepping. 355 356 The residual norms printed by monitoring routines such as 357 SNESDefaultMonitor() (as activated via -snes_monitor) will not be 358 correct, since they are not computed. 359 360 Level: advanced 361 362 .keywords: SNES, nonlinear, line search, cubic 363 364 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 365 SNESSetLineSearch(), SNESNoLineSearch() 366 @*/ 367 int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 368 { 369 int ierr; 370 Scalar mone = -1.0; 371 SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 372 PetscTruth change_y = PETSC_FALSE; 373 374 PetscFunctionBegin; 375 *flag = 0; 376 PLogEventBegin(SNES_LineSearch,snes,x,f,g); 377 ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 378 if (neP->CheckStep) { 379 ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 380 } 381 382 /* don't evaluate function the last time through */ 383 if (snes->iter < snes->max_its-1) { 384 ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 385 } 386 PLogEventEnd(SNES_LineSearch,snes,x,f,g); 387 PetscFunctionReturn(0); 388 } 389 /* -------------------------------------------------------------------------- */ 390 #undef __FUNC__ 391 #define __FUNC__ /*<a name=""></a>*/"SNESCubicLineSearch" 392 /*@C 393 SNESCubicLineSearch - Performs a cubic line search (default line search method). 394 395 Collective on SNES 396 397 Input Parameters: 398 + snes - nonlinear context 399 . lsctx - optional context for line search (not used here) 400 . x - current iterate 401 . f - residual evaluated at x 402 . y - search direction (contains new iterate on output) 403 . w - work vector 404 - fnorm - 2-norm of f 405 406 Output Parameters: 407 + g - residual evaluated at new iterate y 408 . y - new iterate (contains search direction on input) 409 . gnorm - 2-norm of g 410 . ynorm - 2-norm of search length 411 - flag - 0 if line search succeeds; -1 on failure. 412 413 Options Database Key: 414 $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 415 416 Notes: 417 This line search is taken from "Numerical Methods for Unconstrained 418 Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 419 420 Level: advanced 421 422 .keywords: SNES, nonlinear, line search, cubic 423 424 .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 425 @*/ 426 int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 427 { 428 /* 429 Note that for line search purposes we work with with the related 430 minimization problem: 431 min z(x): R^n -> R, 432 where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 433 */ 434 435 PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2; 436 PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 437 #if defined(PETSC_USE_COMPLEX) 438 Scalar cinitslope,clambda; 439 #endif 440 int ierr,count; 441 SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 442 Scalar mone = -1.0,scale; 443 PetscTruth change_y = PETSC_FALSE; 444 445 PetscFunctionBegin; 446 PLogEventBegin(SNES_LineSearch,snes,x,f,g); 447 *flag = 0; 448 alpha = neP->alpha; 449 maxstep = neP->maxstep; 450 steptol = neP->steptol; 451 452 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 453 if (*ynorm < snes->atol) { 454 PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 455 *gnorm = fnorm; 456 ierr = VecCopy(x,y);CHKERRQ(ierr); 457 ierr = VecCopy(f,g);CHKERRQ(ierr); 458 goto theend1; 459 } 460 if (*ynorm > maxstep) { /* Step too big, so scale back */ 461 scale = maxstep/(*ynorm); 462 #if defined(PETSC_USE_COMPLEX) 463 PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm); 464 #else 465 PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm); 466 #endif 467 ierr = VecScale(&scale,y);CHKERRQ(ierr); 468 *ynorm = maxstep; 469 } 470 minlambda = steptol/(*ynorm); 471 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 472 #if defined(PETSC_USE_COMPLEX) 473 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 474 initslope = PetscRealPart(cinitslope); 475 #else 476 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 477 #endif 478 if (initslope > 0.0) initslope = -initslope; 479 if (initslope == 0.0) initslope = -1.0; 480 481 ierr = VecCopy(y,w);CHKERRQ(ierr); 482 ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 483 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 484 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 485 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 486 ierr = VecCopy(w,y);CHKERRQ(ierr); 487 PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 488 goto theend1; 489 } 490 491 /* Fit points with quadratic */ 492 lambda = 1.0; 493 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 494 lambdaprev = lambda; 495 gnormprev = *gnorm; 496 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 497 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 498 else lambda = lambdatemp; 499 ierr = VecCopy(x,w);CHKERRQ(ierr); 500 lambdaneg = -lambda; 501 #if defined(PETSC_USE_COMPLEX) 502 clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 503 #else 504 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 505 #endif 506 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 507 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 508 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 509 ierr = VecCopy(w,y);CHKERRQ(ierr); 510 PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 511 goto theend1; 512 } 513 514 /* Fit points with cubic */ 515 count = 1; 516 while (1) { 517 if (lambda <= minlambda) { /* bad luck; use full step */ 518 PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 519 PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 520 ierr = VecCopy(w,y);CHKERRQ(ierr); 521 *flag = -1; break; 522 } 523 t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 524 t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 525 a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 526 b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 527 d = b*b - 3*a*initslope; 528 if (d < 0.0) d = 0.0; 529 if (a == 0.0) { 530 lambdatemp = -initslope/(2.0*b); 531 } else { 532 lambdatemp = (-b + sqrt(d))/(3.0*a); 533 } 534 lambdaprev = lambda; 535 gnormprev = *gnorm; 536 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 537 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 538 else lambda = lambdatemp; 539 ierr = VecCopy(x,w);CHKERRQ(ierr); 540 lambdaneg = -lambda; 541 #if defined(PETSC_USE_COMPLEX) 542 clambda = lambdaneg; 543 ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 544 #else 545 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 546 #endif 547 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 548 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 549 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 550 ierr = VecCopy(w,y);CHKERRQ(ierr); 551 PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 552 goto theend1; 553 } else { 554 PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 555 } 556 count++; 557 } 558 theend1: 559 /* Optional user-defined check for line search step validity */ 560 if (neP->CheckStep) { 561 ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 562 if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 563 ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 564 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 565 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 566 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 567 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 568 } 569 } 570 PLogEventEnd(SNES_LineSearch,snes,x,f,g); 571 PetscFunctionReturn(0); 572 } 573 /* -------------------------------------------------------------------------- */ 574 #undef __FUNC__ 575 #define __FUNC__ /*<a name=""></a>*/"SNESQuadraticLineSearch" 576 /*@C 577 SNESQuadraticLineSearch - Performs a quadratic line search. 578 579 Collective on SNES and Vec 580 581 Input Parameters: 582 + snes - the SNES context 583 . lsctx - optional context for line search (not used here) 584 . x - current iterate 585 . f - residual evaluated at x 586 . y - search direction (contains new iterate on output) 587 . w - work vector 588 - fnorm - 2-norm of f 589 590 Output Parameters: 591 + g - residual evaluated at new iterate y 592 . y - new iterate (contains search direction on input) 593 . gnorm - 2-norm of g 594 . ynorm - 2-norm of search length 595 - flag - 0 if line search succeeds; -1 on failure. 596 597 Options Database Key: 598 . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 599 600 Notes: 601 Use SNESSetLineSearch() to set this routine within the SNESEQLS method. 602 603 Level: advanced 604 605 .keywords: SNES, nonlinear, quadratic, line search 606 607 .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 608 @*/ 609 int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 610 { 611 /* 612 Note that for line search purposes we work with with the related 613 minimization problem: 614 min z(x): R^n -> R, 615 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 616 */ 617 PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 618 #if defined(PETSC_USE_COMPLEX) 619 Scalar cinitslope,clambda; 620 #endif 621 int ierr,count; 622 SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 623 Scalar mone = -1.0,scale; 624 PetscTruth change_y = PETSC_FALSE; 625 626 PetscFunctionBegin; 627 PLogEventBegin(SNES_LineSearch,snes,x,f,g); 628 *flag = 0; 629 alpha = neP->alpha; 630 maxstep = neP->maxstep; 631 steptol = neP->steptol; 632 633 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 634 if (*ynorm < snes->atol) { 635 PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 636 *gnorm = fnorm; 637 ierr = VecCopy(x,y);CHKERRQ(ierr); 638 ierr = VecCopy(f,g);CHKERRQ(ierr); 639 goto theend2; 640 } 641 if (*ynorm > maxstep) { /* Step too big, so scale back */ 642 scale = maxstep/(*ynorm); 643 ierr = VecScale(&scale,y);CHKERRQ(ierr); 644 *ynorm = maxstep; 645 } 646 minlambda = steptol/(*ynorm); 647 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 648 #if defined(PETSC_USE_COMPLEX) 649 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 650 initslope = PetscRealPart(cinitslope); 651 #else 652 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 653 #endif 654 if (initslope > 0.0) initslope = -initslope; 655 if (initslope == 0.0) initslope = -1.0; 656 657 ierr = VecCopy(y,w);CHKERRQ(ierr); 658 ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 659 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 660 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 661 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 662 ierr = VecCopy(w,y);CHKERRQ(ierr); 663 PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 664 goto theend2; 665 } 666 667 /* Fit points with quadratic */ 668 lambda = 1.0; 669 count = 1; 670 while (1) { 671 if (lambda <= minlambda) { /* bad luck; use full step */ 672 PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 673 PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 674 ierr = VecCopy(w,y);CHKERRQ(ierr); 675 *flag = -1; break; 676 } 677 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 678 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 679 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 680 else lambda = lambdatemp; 681 ierr = VecCopy(x,w);CHKERRQ(ierr); 682 lambdaneg = -lambda; 683 #if defined(PETSC_USE_COMPLEX) 684 clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 685 #else 686 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 687 #endif 688 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 689 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 690 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 691 ierr = VecCopy(w,y);CHKERRQ(ierr); 692 PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 693 break; 694 } 695 count++; 696 } 697 theend2: 698 /* Optional user-defined check for line search step validity */ 699 if (neP->CheckStep) { 700 ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 701 if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 702 ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 703 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 704 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 705 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 706 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 707 } 708 } 709 PLogEventEnd(SNES_LineSearch,snes,x,f,g); 710 PetscFunctionReturn(0); 711 } 712 /* -------------------------------------------------------------------------- */ 713 #undef __FUNC__ 714 #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch" 715 /*@C 716 SNESSetLineSearch - Sets the line search routine to be used 717 by the method SNESEQLS. 718 719 Input Parameters: 720 + snes - nonlinear context obtained from SNESCreate() 721 . lsctx - optional user-defined context for use by line search 722 - func - pointer to int function 723 724 Collective on SNES 725 726 Available Routines: 727 + SNESCubicLineSearch() - default line search 728 . SNESQuadraticLineSearch() - quadratic line search 729 . SNESNoLineSearch() - the full Newton step (actually not a line search) 730 - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 731 732 Options Database Keys: 733 + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 734 . -snes_eq_ls_alpha <alpha> - Sets alpha 735 . -snes_eq_ls_maxstep <max> - Sets maxstep 736 - -snes_eq_ls_steptol <steptol> - Sets steptol 737 738 Calling sequence of func: 739 .vb 740 func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 741 PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag) 742 .ve 743 744 Input parameters for func: 745 + snes - nonlinear context 746 . lsctx - optional user-defined context for line search 747 . x - current iterate 748 . f - residual evaluated at x 749 . y - search direction (contains new iterate on output) 750 . w - work vector 751 - fnorm - 2-norm of f 752 753 Output parameters for func: 754 + g - residual evaluated at new iterate y 755 . y - new iterate (contains search direction on input) 756 . gnorm - 2-norm of g 757 . ynorm - 2-norm of search length 758 - flag - set to 0 if the line search succeeds; a nonzero integer 759 on failure. 760 761 Level: advanced 762 763 .keywords: SNES, nonlinear, set, line search, routine 764 765 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), 766 SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 767 @*/ 768 int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 769 { 770 int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*); 771 772 PetscFunctionBegin; 773 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 774 if (f) { 775 ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 776 } 777 PetscFunctionReturn(0); 778 } 779 /* -------------------------------------------------------------------------- */ 780 EXTERN_C_BEGIN 781 #undef __FUNC__ 782 #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch_LS" 783 int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 784 double,double*,double*,int*),void *lsctx) 785 { 786 PetscFunctionBegin; 787 ((SNES_EQ_LS *)(snes->data))->LineSearch = func; 788 ((SNES_EQ_LS *)(snes->data))->lsP = lsctx; 789 PetscFunctionReturn(0); 790 } 791 EXTERN_C_END 792 /* -------------------------------------------------------------------------- */ 793 #undef __FUNC__ 794 #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck" 795 /*@C 796 SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 797 by the line search routine in the Newton-based method SNESEQLS. 798 799 Input Parameters: 800 + snes - nonlinear context obtained from SNESCreate() 801 . func - pointer to int function 802 - checkctx - optional user-defined context for use by step checking routine 803 804 Collective on SNES 805 806 Calling sequence of func: 807 .vb 808 int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 809 .ve 810 where func returns an error code of 0 on success and a nonzero 811 on failure. 812 813 Input parameters for func: 814 + snes - nonlinear context 815 . checkctx - optional user-defined context for use by step checking routine 816 - x - current candidate iterate 817 818 Output parameters for func: 819 + x - current iterate (possibly modified) 820 - flag - flag indicating whether x has been modified (either 821 PETSC_TRUE of PETSC_FALSE) 822 823 Level: advanced 824 825 Notes: 826 SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 827 iterate computed by the line search checking routine. In particular, 828 these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 829 to the checking routine, and then (3) compute the corresponding nonlinear 830 function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 831 832 SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 833 new iterate computed by the line search checking routine. In particular, 834 these routines (1) compute a candidate iterate u_{i+1} as well as a 835 candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 836 routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 837 were made to the candidate iterate in the checking routine (as indicated 838 by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 839 very costly, so use this feature with caution! 840 841 .keywords: SNES, nonlinear, set, line search check, step check, routine 842 843 .seealso: SNESSetLineSearch() 844 @*/ 845 int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 846 { 847 int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 848 849 PetscFunctionBegin; 850 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr); 851 if (f) { 852 ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 853 } 854 PetscFunctionReturn(0); 855 } 856 /* -------------------------------------------------------------------------- */ 857 EXTERN_C_BEGIN 858 #undef __FUNC__ 859 #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck_LS" 860 int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 861 { 862 PetscFunctionBegin; 863 ((SNES_EQ_LS *)(snes->data))->CheckStep = func; 864 ((SNES_EQ_LS *)(snes->data))->checkP = checkctx; 865 PetscFunctionReturn(0); 866 } 867 EXTERN_C_END 868 /* -------------------------------------------------------------------------- */ 869 /* 870 SNESPrintHelp_EQ_LS - Prints all options for the SNESEQLS method. 871 872 Input Parameter: 873 . snes - the SNES context 874 875 Application Interface Routine: SNESPrintHelp() 876 */ 877 #undef __FUNC__ 878 #define __FUNC__ /*<a name=""></a>*/"SNESPrintHelp_EQ_LS" 879 static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 880 { 881 SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 882 int ierr; 883 884 PetscFunctionBegin; 885 ierr = (*PetscHelpPrintf)(snes->comm," method SNESEQLS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr); 886 ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr); 887 ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr); 888 ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr); 889 ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr); 890 PetscFunctionReturn(0); 891 } 892 /* -------------------------------------------------------------------------- */ 893 /* 894 SNESView_EQ_LS - Prints info from the SNESEQLS data structure. 895 896 Input Parameters: 897 . SNES - the SNES context 898 . viewer - visualization context 899 900 Application Interface Routine: SNESView() 901 */ 902 #undef __FUNC__ 903 #define __FUNC__ /*<a name=""></a>*/"SNESView_EQ_LS" 904 static int SNESView_EQ_LS(SNES snes,Viewer viewer) 905 { 906 SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 907 char *cstr; 908 int ierr; 909 PetscTruth isascii; 910 911 PetscFunctionBegin; 912 ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 913 if (isascii) { 914 if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 915 else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 916 else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 917 else cstr = "unknown"; 918 ierr = ViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 919 ierr = ViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 920 } else { 921 SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 922 } 923 PetscFunctionReturn(0); 924 } 925 /* -------------------------------------------------------------------------- */ 926 /* 927 SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method. 928 929 Input Parameter: 930 . snes - the SNES context 931 932 Application Interface Routine: SNESSetFromOptions() 933 */ 934 #undef __FUNC__ 935 #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions_EQ_LS" 936 static int SNESSetFromOptions_EQ_LS(SNES snes) 937 { 938 SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 939 char ver[16]; 940 double tmp; 941 int ierr; 942 PetscTruth flg; 943 944 PetscFunctionBegin; 945 ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp,&flg);CHKERRQ(ierr); 946 if (flg) { 947 ls->alpha = tmp; 948 } 949 ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp,&flg);CHKERRQ(ierr); 950 if (flg) { 951 ls->maxstep = tmp; 952 } 953 ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp,&flg);CHKERRQ(ierr); 954 if (flg) { 955 ls->steptol = tmp; 956 } 957 ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16,&flg);CHKERRQ(ierr); 958 if (flg) { 959 PetscTruth isbasic,isnonorms,isquad,iscubic; 960 961 ierr = PetscStrcmp(ver,"basic",&isbasic);CHKERRQ(ierr); 962 ierr = PetscStrcmp(ver,"basicnonorms",&isnonorms);CHKERRQ(ierr); 963 ierr = PetscStrcmp(ver,"quadratic",&isquad);CHKERRQ(ierr); 964 ierr = PetscStrcmp(ver,"cubic",&iscubic);CHKERRQ(ierr); 965 966 if (isbasic) { 967 ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 968 } else if (isnonorms) { 969 ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 970 } else if (isquad) { 971 ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 972 } else if (iscubic) { 973 ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 974 } 975 else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 976 } 977 PetscFunctionReturn(0); 978 } 979 /* -------------------------------------------------------------------------- */ 980 /* 981 SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method, 982 SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver 983 context, SNES, that was created within SNESCreate(). 984 985 986 Input Parameter: 987 . snes - the SNES context 988 989 Application Interface Routine: SNESCreate() 990 */ 991 EXTERN_C_BEGIN 992 #undef __FUNC__ 993 #define __FUNC__ /*<a name=""></a>*/"SNESCreate_EQ_LS" 994 int SNESCreate_EQ_LS(SNES snes) 995 { 996 int ierr; 997 SNES_EQ_LS *neP; 998 999 PetscFunctionBegin; 1000 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1001 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 1002 } 1003 1004 snes->setup = SNESSetUp_EQ_LS; 1005 snes->solve = SNESSolve_EQ_LS; 1006 snes->destroy = SNESDestroy_EQ_LS; 1007 snes->converged = SNESConverged_EQ_LS; 1008 snes->printhelp = SNESPrintHelp_EQ_LS; 1009 snes->setfromoptions = SNESSetFromOptions_EQ_LS; 1010 snes->view = SNESView_EQ_LS; 1011 snes->nwork = 0; 1012 1013 neP = PetscNew(SNES_EQ_LS);CHKPTRQ(neP); 1014 PLogObjectMemory(snes,sizeof(SNES_EQ_LS)); 1015 snes->data = (void*)neP; 1016 neP->alpha = 1.e-4; 1017 neP->maxstep = 1.e8; 1018 neP->steptol = 1.e-12; 1019 neP->LineSearch = SNESCubicLineSearch; 1020 neP->lsP = PETSC_NULL; 1021 neP->CheckStep = PETSC_NULL; 1022 neP->checkP = PETSC_NULL; 1023 1024 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr); 1025 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 1026 1027 PetscFunctionReturn(0); 1028 } 1029 EXTERN_C_END 1030 1031 1032 1033