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