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