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