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