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