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