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