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