1 /*$Id: ls.c,v 1.149 1999/12/20 22:57:33 bsmith Exp bsmith $*/ 2 3 #include "src/snes/impls/ls/ls.h" 4 5 /* -------------------------------------------------------------------- 6 7 This file implements a truncated Newton method with a line search, 8 for solving a system of nonlinear equations, using the SLES, Vec, 9 and Mat interfaces for linear solvers, vectors, and matrices, 10 respectively. 11 12 The following basic routines are required for each nonlinear solver: 13 SNESCreate_XXX() - Creates a nonlinear solver context 14 SNESSetFromOptions_XXX() - Sets runtime options 15 SNESSolve_XXX() - Solves the nonlinear system 16 SNESDestroy_XXX() - Destroys the nonlinear solver context 17 The suffix "_XXX" denotes a particular implementation, in this case 18 we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving 19 systems of nonlinear equations (EQ) with a line search (LS) method. 20 These routines are actually called via the common user interface 21 routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 22 SNESDestroy(), so the application code interface remains identical 23 for all nonlinear solvers. 24 25 Another key routine is: 26 SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 27 by setting data structures and options. The interface routine SNESSetUp() 28 is not usually called directly by the user, but instead is called by 29 SNESSolve() if necessary. 30 31 Additional basic routines are: 32 SNESPrintHelp_XXX() - Prints nonlinear solver runtime options 33 SNESView_XXX() - Prints details of runtime options that 34 have actually been used. 35 These are called by application codes via the interface routines 36 SNESPrintHelp() and SNESView(). 37 38 The various types of solvers (preconditioners, Krylov subspace methods, 39 nonlinear solvers, timesteppers) are all organized similarly, so the 40 above description applies to these categories also. 41 42 -------------------------------------------------------------------- */ 43 /* 44 SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton 45 method with a line search. 46 47 Input Parameters: 48 . snes - the SNES context 49 50 Output Parameter: 51 . outits - number of iterations until termination 52 53 Application Interface Routine: SNESSolve() 54 55 Notes: 56 This implements essentially a truncated Newton method with a 57 line search. By default a cubic backtracking line search 58 is employed, as described in the text "Numerical Methods for 59 Unconstrained Optimization and Nonlinear Equations" by Dennis 60 and Schnabel. 61 */ 62 #undef __FUNC__ 63 #define __FUNC__ "SNESSolve_EQ_LS" 64 int SNESSolve_EQ_LS(SNES snes,int *outits) 65 { 66 SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 67 int maxits,i,ierr,lits,lsfail; 68 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 69 PetscReal fnorm,gnorm,xnorm,ynorm; 70 Vec Y,X,F,G,W,TMP; 71 72 PetscFunctionBegin; 73 snes->reason = SNES_CONVERGED_ITERATING; 74 75 maxits = snes->max_its; /* maximum number of iterations */ 76 X = snes->vec_sol; /* solution vector */ 77 F = snes->vec_func; /* residual vector */ 78 Y = snes->work[0]; /* work vectors */ 79 G = snes->work[1]; 80 W = snes->work[2]; 81 82 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 83 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 84 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 85 snes->iter = 0; 86 snes->norm = fnorm; 87 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 88 SNESLogConvHistory(snes,fnorm,0); 89 SNESMonitor(snes,0,fnorm); 90 91 if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; 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 += 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) { 113 snes->nfailures++; 114 snes->reason = SNES_DIVERGED_LS_FAILURE; 115 break; 116 } 117 118 TMP = F; F = G; snes->vec_func_always = F; G = TMP; 119 TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP; 120 fnorm = gnorm; 121 122 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 123 snes->iter = i+1; 124 snes->norm = fnorm; 125 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 126 SNESLogConvHistory(snes,fnorm,lits); 127 SNESMonitor(snes,i+1,fnorm); 128 129 /* Test for convergence */ 130 if (snes->converged) { 131 ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 132 ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 133 if (snes->reason) { 134 break; 135 } 136 } 137 } 138 if (X != snes->vec_sol) { 139 ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 140 snes->vec_sol_always = snes->vec_sol; 141 snes->vec_func_always = snes->vec_func; 142 } 143 if (i == maxits) { 144 PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits); 145 i--; 146 snes->reason = SNES_DIVERGED_MAX_IT; 147 } 148 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 149 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 150 *outits = i+1; 151 PetscFunctionReturn(0); 152 } 153 /* -------------------------------------------------------------------------- */ 154 /* 155 SNESSetUp_EQ_LS - Sets up the internal data structures for the later use 156 of the SNESEQLS nonlinear solver. 157 158 Input Parameter: 159 . snes - the SNES context 160 . x - the solution vector 161 162 Application Interface Routine: SNESSetUp() 163 164 Notes: 165 For basic use of the SNES solvers the user need not explicitly call 166 SNESSetUp(), since these actions will automatically occur during 167 the call to SNESSolve(). 168 */ 169 #undef __FUNC__ 170 #define __FUNC__ "SNESSetUp_EQ_LS" 171 int SNESSetUp_EQ_LS(SNES snes) 172 { 173 int ierr; 174 175 PetscFunctionBegin; 176 snes->nwork = 4; 177 ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 178 PLogObjectParents(snes,snes->nwork,snes->work); 179 snes->vec_sol_update_always = snes->work[3]; 180 PetscFunctionReturn(0); 181 } 182 /* -------------------------------------------------------------------------- */ 183 /* 184 SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created 185 with SNESCreate_EQ_LS(). 186 187 Input Parameter: 188 . snes - the SNES context 189 190 Application Interface Routine: SNESDestroy() 191 */ 192 #undef __FUNC__ 193 #define __FUNC__ "SNESDestroy_EQ_LS" 194 int SNESDestroy_EQ_LS(SNES snes) 195 { 196 int ierr; 197 198 PetscFunctionBegin; 199 if (snes->nwork) { 200 ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 201 } 202 ierr = PetscFree(snes->data);CHKERRQ(ierr); 203 PetscFunctionReturn(0); 204 } 205 /* -------------------------------------------------------------------------- */ 206 #undef __FUNC__ 207 #define __FUNC__ "SNESNoLineSearch" 208 209 /*@C 210 SNESNoLineSearch - This routine is not a line search at all; 211 it simply uses the full Newton step. Thus, this routine is intended 212 to serve as a template and is not recommended for general use. 213 214 Collective on SNES and Vec 215 216 Input Parameters: 217 + snes - nonlinear context 218 . lsctx - optional context for line search (not used here) 219 . x - current iterate 220 . f - residual evaluated at x 221 . y - search direction (contains new iterate on output) 222 . w - work vector 223 - fnorm - 2-norm of f 224 225 Output Parameters: 226 + g - residual evaluated at new iterate y 227 . y - new iterate (contains search direction on input) 228 . gnorm - 2-norm of g 229 . ynorm - 2-norm of search length 230 - flag - set to 0, indicating a successful line search 231 232 Options Database Key: 233 . -snes_eq_ls basic - Activates SNESNoLineSearch() 234 235 Level: advanced 236 237 .keywords: SNES, nonlinear, line search, cubic 238 239 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 240 SNESSetLineSearch(), SNESNoLineSearchNoNorms() 241 @*/ 242 int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 243 PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 244 { 245 int ierr; 246 Scalar mone = -1.0; 247 SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 248 PetscTruth change_y = PETSC_FALSE; 249 250 PetscFunctionBegin; 251 *flag = 0; 252 PLogEventBegin(SNES_LineSearch,snes,x,f,g); 253 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 254 ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 255 if (neP->CheckStep) { 256 ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 257 } 258 ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 259 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 260 PLogEventEnd(SNES_LineSearch,snes,x,f,g); 261 PetscFunctionReturn(0); 262 } 263 /* -------------------------------------------------------------------------- */ 264 #undef __FUNC__ 265 #define __FUNC__ "SNESNoLineSearchNoNorms" 266 267 /*@C 268 SNESNoLineSearchNoNorms - This routine is not a line search at 269 all; it simply uses the full Newton step. This version does not 270 even compute the norm of the function or search direction; this 271 is intended only when you know the full step is fine and are 272 not checking for convergence of the nonlinear iteration (for 273 example, you are running always for a fixed number of Newton steps). 274 275 Collective on SNES and Vec 276 277 Input Parameters: 278 + snes - nonlinear context 279 . lsctx - optional context for line search (not used here) 280 . x - current iterate 281 . f - residual evaluated at x 282 . y - search direction (contains new iterate on output) 283 . w - work vector 284 - fnorm - 2-norm of f 285 286 Output Parameters: 287 + g - residual evaluated at new iterate y 288 . gnorm - not changed 289 . ynorm - not changed 290 - flag - set to 0, indicating a successful line search 291 292 Options Database Key: 293 . -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 294 295 Notes: 296 SNESNoLineSearchNoNorms() must be used in conjunction with 297 either the options 298 $ -snes_no_convergence_test -snes_max_it <its> 299 or alternatively a user-defined custom test set via 300 SNESSetConvergenceTest(); or a -snes_max_it of 1, 301 otherwise, the SNES solver will generate an error. 302 303 During the final iteration this will not evaluate the function at 304 the solution point. This is to save a function evaluation while 305 using pseudo-timestepping. 306 307 The residual norms printed by monitoring routines such as 308 SNESDefaultMonitor() (as activated via -snes_monitor) will not be 309 correct, since they are not computed. 310 311 Level: advanced 312 313 .keywords: SNES, nonlinear, line search, cubic 314 315 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 316 SNESSetLineSearch(), SNESNoLineSearch() 317 @*/ 318 int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 319 PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 320 { 321 int ierr; 322 Scalar mone = -1.0; 323 SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 324 PetscTruth change_y = PETSC_FALSE; 325 326 PetscFunctionBegin; 327 *flag = 0; 328 PLogEventBegin(SNES_LineSearch,snes,x,f,g); 329 ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr); /* y <- y - x */ 330 if (neP->CheckStep) { 331 ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 332 } 333 334 /* don't evaluate function the last time through */ 335 if (snes->iter < snes->max_its-1) { 336 ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y) */ 337 } 338 PLogEventEnd(SNES_LineSearch,snes,x,f,g); 339 PetscFunctionReturn(0); 340 } 341 /* -------------------------------------------------------------------------- */ 342 #undef __FUNC__ 343 #define __FUNC__ "SNESCubicLineSearch" 344 /*@C 345 SNESCubicLineSearch - Performs a cubic line search (default line search method). 346 347 Collective on SNES 348 349 Input Parameters: 350 + snes - nonlinear context 351 . lsctx - optional context for line search (not used here) 352 . x - current iterate 353 . f - residual evaluated at x 354 . y - search direction (contains new iterate on output) 355 . w - work vector 356 - fnorm - 2-norm of f 357 358 Output Parameters: 359 + g - residual evaluated at new iterate y 360 . y - new iterate (contains search direction on input) 361 . gnorm - 2-norm of g 362 . ynorm - 2-norm of search length 363 - flag - 0 if line search succeeds; -1 on failure. 364 365 Options Database Key: 366 $ -snes_eq_ls cubic - Activates SNESCubicLineSearch() 367 368 Notes: 369 This line search is taken from "Numerical Methods for Unconstrained 370 Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 371 372 Level: advanced 373 374 .keywords: SNES, nonlinear, line search, cubic 375 376 .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 377 @*/ 378 int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 379 PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 380 { 381 /* 382 Note that for line search purposes we work with with the related 383 minimization problem: 384 min z(x): R^n -> R, 385 where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 386 */ 387 388 PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2; 389 PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 390 #if defined(PETSC_USE_COMPLEX) 391 Scalar cinitslope,clambda; 392 #endif 393 int ierr,count; 394 SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 395 Scalar mone = -1.0,scale; 396 PetscTruth change_y = PETSC_FALSE; 397 398 PetscFunctionBegin; 399 PLogEventBegin(SNES_LineSearch,snes,x,f,g); 400 *flag = 0; 401 alpha = neP->alpha; 402 maxstep = neP->maxstep; 403 steptol = neP->steptol; 404 405 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 406 if (*ynorm < snes->atol) { 407 PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n"); 408 *gnorm = fnorm; 409 ierr = VecCopy(x,y);CHKERRQ(ierr); 410 ierr = VecCopy(f,g);CHKERRQ(ierr); 411 goto theend1; 412 } 413 if (*ynorm > maxstep) { /* Step too big, so scale back */ 414 scale = maxstep/(*ynorm); 415 #if defined(PETSC_USE_COMPLEX) 416 PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscRealPart(scale)); 417 #else 418 PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale); 419 #endif 420 ierr = VecScale(&scale,y);CHKERRQ(ierr); 421 *ynorm = maxstep; 422 } 423 minlambda = steptol/(*ynorm); 424 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 425 #if defined(PETSC_USE_COMPLEX) 426 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 427 initslope = PetscRealPart(cinitslope); 428 #else 429 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 430 #endif 431 if (initslope > 0.0) initslope = -initslope; 432 if (initslope == 0.0) initslope = -1.0; 433 434 ierr = VecCopy(y,w);CHKERRQ(ierr); 435 ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 436 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 437 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 438 if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 439 ierr = VecCopy(w,y);CHKERRQ(ierr); 440 PLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 441 goto theend1; 442 } 443 444 /* Fit points with quadratic */ 445 lambda = 1.0; 446 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 447 lambdaprev = lambda; 448 gnormprev = *gnorm; 449 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 450 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 451 else lambda = lambdatemp; 452 ierr = VecCopy(x,w);CHKERRQ(ierr); 453 lambdaneg = -lambda; 454 #if defined(PETSC_USE_COMPLEX) 455 clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 456 #else 457 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 458 #endif 459 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 460 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 461 if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 462 ierr = VecCopy(w,y);CHKERRQ(ierr); 463 PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda); 464 goto theend1; 465 } 466 467 /* Fit points with cubic */ 468 count = 1; 469 while (1) { 470 if (lambda <= minlambda) { /* bad luck; use full step */ 471 PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 472 PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 473 fnorm,*gnorm,*ynorm,lambda,initslope); 474 ierr = VecCopy(w,y);CHKERRQ(ierr); 475 *flag = -1; break; 476 } 477 t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope; 478 t2 = (gnormprev*gnormprev - fnorm*fnorm)*0.5 - lambdaprev*initslope; 479 a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 480 b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 481 d = b*b - 3*a*initslope; 482 if (d < 0.0) d = 0.0; 483 if (a == 0.0) { 484 lambdatemp = -initslope/(2.0*b); 485 } else { 486 lambdatemp = (-b + sqrt(d))/(3.0*a); 487 } 488 lambdaprev = lambda; 489 gnormprev = *gnorm; 490 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 491 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 492 else lambda = lambdatemp; 493 ierr = VecCopy(x,w);CHKERRQ(ierr); 494 lambdaneg = -lambda; 495 #if defined(PETSC_USE_COMPLEX) 496 clambda = lambdaneg; 497 ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 498 #else 499 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 500 #endif 501 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 502 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 503 if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */ 504 ierr = VecCopy(w,y);CHKERRQ(ierr); 505 PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda); 506 goto theend1; 507 } else { 508 PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%g\n",lambda); 509 } 510 count++; 511 } 512 theend1: 513 /* Optional user-defined check for line search step validity */ 514 if (neP->CheckStep) { 515 ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 516 if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 517 ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 518 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 519 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 520 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 521 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 522 } 523 } 524 PLogEventEnd(SNES_LineSearch,snes,x,f,g); 525 PetscFunctionReturn(0); 526 } 527 /* -------------------------------------------------------------------------- */ 528 #undef __FUNC__ 529 #define __FUNC__ "SNESQuadraticLineSearch" 530 /*@C 531 SNESQuadraticLineSearch - Performs a quadratic line search. 532 533 Collective on SNES and Vec 534 535 Input Parameters: 536 + snes - the SNES context 537 . lsctx - optional context for line search (not used here) 538 . x - current iterate 539 . f - residual evaluated at x 540 . y - search direction (contains new iterate on output) 541 . w - work vector 542 - fnorm - 2-norm of f 543 544 Output Parameters: 545 + g - residual evaluated at new iterate y 546 . y - new iterate (contains search direction on input) 547 . gnorm - 2-norm of g 548 . ynorm - 2-norm of search length 549 - flag - 0 if line search succeeds; -1 on failure. 550 551 Options Database Key: 552 . -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch() 553 554 Notes: 555 Use SNESSetLineSearch() to set this routine within the SNESEQLS method. 556 557 Level: advanced 558 559 .keywords: SNES, nonlinear, quadratic, line search 560 561 .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 562 @*/ 563 int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 564 PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 565 { 566 /* 567 Note that for line search purposes we work with with the related 568 minimization problem: 569 min z(x): R^n -> R, 570 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 571 */ 572 PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 573 #if defined(PETSC_USE_COMPLEX) 574 Scalar cinitslope,clambda; 575 #endif 576 int ierr,count; 577 SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data; 578 Scalar mone = -1.0,scale; 579 PetscTruth change_y = PETSC_FALSE; 580 581 PetscFunctionBegin; 582 PLogEventBegin(SNES_LineSearch,snes,x,f,g); 583 *flag = 0; 584 alpha = neP->alpha; 585 maxstep = neP->maxstep; 586 steptol = neP->steptol; 587 588 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 589 if (*ynorm < snes->atol) { 590 PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 591 *gnorm = fnorm; 592 ierr = VecCopy(x,y);CHKERRQ(ierr); 593 ierr = VecCopy(f,g);CHKERRQ(ierr); 594 goto theend2; 595 } 596 if (*ynorm > maxstep) { /* Step too big, so scale back */ 597 scale = maxstep/(*ynorm); 598 ierr = VecScale(&scale,y);CHKERRQ(ierr); 599 *ynorm = maxstep; 600 } 601 minlambda = steptol/(*ynorm); 602 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 603 #if defined(PETSC_USE_COMPLEX) 604 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 605 initslope = PetscRealPart(cinitslope); 606 #else 607 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 608 #endif 609 if (initslope > 0.0) initslope = -initslope; 610 if (initslope == 0.0) initslope = -1.0; 611 612 ierr = VecCopy(y,w);CHKERRQ(ierr); 613 ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 614 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 615 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 616 if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */ 617 ierr = VecCopy(w,y);CHKERRQ(ierr); 618 PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 619 goto theend2; 620 } 621 622 /* Fit points with quadratic */ 623 lambda = 1.0; 624 count = 1; 625 while (1) { 626 if (lambda <= minlambda) { /* bad luck; use full step */ 627 PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 628 PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n", 629 fnorm,*gnorm,*ynorm,lambda,initslope); 630 ierr = VecCopy(w,y);CHKERRQ(ierr); 631 *flag = -1; break; 632 } 633 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 634 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 635 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 636 else lambda = lambdatemp; 637 ierr = VecCopy(x,w);CHKERRQ(ierr); 638 lambdaneg = -lambda; 639 #if defined(PETSC_USE_COMPLEX) 640 clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 641 #else 642 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 643 #endif 644 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 645 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 646 if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */ 647 ierr = VecCopy(w,y);CHKERRQ(ierr); 648 PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 649 break; 650 } 651 count++; 652 } 653 theend2: 654 /* Optional user-defined check for line search step validity */ 655 if (neP->CheckStep) { 656 ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 657 if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 658 ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 659 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 660 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 661 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 662 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 663 } 664 } 665 PLogEventEnd(SNES_LineSearch,snes,x,f,g); 666 PetscFunctionReturn(0); 667 } 668 /* -------------------------------------------------------------------------- */ 669 #undef __FUNC__ 670 #define __FUNC__ "SNESSetLineSearch" 671 /*@C 672 SNESSetLineSearch - Sets the line search routine to be used 673 by the method SNESEQLS. 674 675 Input Parameters: 676 + snes - nonlinear context obtained from SNESCreate() 677 . lsctx - optional user-defined context for use by line search 678 - func - pointer to int function 679 680 Collective on SNES 681 682 Available Routines: 683 + SNESCubicLineSearch() - default line search 684 . SNESQuadraticLineSearch() - quadratic line search 685 . SNESNoLineSearch() - the full Newton step (actually not a line search) 686 - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 687 688 Options Database Keys: 689 + -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 690 . -snes_eq_ls_alpha <alpha> - Sets alpha 691 . -snes_eq_ls_maxstep <max> - Sets maxstep 692 - -snes_eq_ls_steptol <steptol> - Sets steptol 693 694 Calling sequence of func: 695 .vb 696 func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 697 PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag) 698 .ve 699 700 Input parameters for func: 701 + snes - nonlinear context 702 . lsctx - optional user-defined context for line search 703 . x - current iterate 704 . f - residual evaluated at x 705 . y - search direction (contains new iterate on output) 706 . w - work vector 707 - fnorm - 2-norm of f 708 709 Output parameters for func: 710 + g - residual evaluated at new iterate y 711 . y - new iterate (contains search direction on input) 712 . gnorm - 2-norm of g 713 . ynorm - 2-norm of search length 714 - flag - set to 0 if the line search succeeds; a nonzero integer 715 on failure. 716 717 Level: advanced 718 719 .keywords: SNES, nonlinear, set, line search, routine 720 721 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 722 @*/ 723 int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 724 { 725 int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*); 726 727 PetscFunctionBegin; 728 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr); 729 if (f) { 730 ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 731 } 732 PetscFunctionReturn(0); 733 } 734 /* -------------------------------------------------------------------------- */ 735 EXTERN_C_BEGIN 736 #undef __FUNC__ 737 #define __FUNC__ "SNESSetLineSearch_LS" 738 int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec, 739 double,double*,double*,int*),void *lsctx) 740 { 741 PetscFunctionBegin; 742 ((SNES_EQ_LS *)(snes->data))->LineSearch = func; 743 ((SNES_EQ_LS *)(snes->data))->lsP = lsctx; 744 PetscFunctionReturn(0); 745 } 746 EXTERN_C_END 747 /* -------------------------------------------------------------------------- */ 748 #undef __FUNC__ 749 #define __FUNC__ "SNESSetLineSearchCheck" 750 /*@C 751 SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 752 by the line search routine in the Newton-based method SNESEQLS. 753 754 Input Parameters: 755 + snes - nonlinear context obtained from SNESCreate() 756 . func - pointer to int function 757 - checkctx - optional user-defined context for use by step checking routine 758 759 Collective on SNES 760 761 Calling sequence of func: 762 .vb 763 int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 764 .ve 765 where func returns an error code of 0 on success and a nonzero 766 on failure. 767 768 Input parameters for func: 769 + snes - nonlinear context 770 . checkctx - optional user-defined context for use by step checking routine 771 - x - current candidate iterate 772 773 Output parameters for func: 774 + x - current iterate (possibly modified) 775 - flag - flag indicating whether x has been modified (either 776 PETSC_TRUE of PETSC_FALSE) 777 778 Level: advanced 779 780 Notes: 781 SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 782 iterate computed by the line search checking routine. In particular, 783 these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 784 to the checking routine, and then (3) compute the corresponding nonlinear 785 function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 786 787 SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 788 new iterate computed by the line search checking routine. In particular, 789 these routines (1) compute a candidate iterate u_{i+1} as well as a 790 candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 791 routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 792 were made to the candidate iterate in the checking routine (as indicated 793 by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 794 very costly, so use this feature with caution! 795 796 .keywords: SNES, nonlinear, set, line search check, step check, routine 797 798 .seealso: SNESSetLineSearch() 799 @*/ 800 int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 801 { 802 int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 803 804 PetscFunctionBegin; 805 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr); 806 if (f) { 807 ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 808 } 809 PetscFunctionReturn(0); 810 } 811 /* -------------------------------------------------------------------------- */ 812 EXTERN_C_BEGIN 813 #undef __FUNC__ 814 #define __FUNC__ "SNESSetLineSearchCheck_LS" 815 int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 816 { 817 PetscFunctionBegin; 818 ((SNES_EQ_LS *)(snes->data))->CheckStep = func; 819 ((SNES_EQ_LS *)(snes->data))->checkP = checkctx; 820 PetscFunctionReturn(0); 821 } 822 EXTERN_C_END 823 /* -------------------------------------------------------------------------- */ 824 /* 825 SNESPrintHelp_EQ_LS - Prints all options for the SNESEQLS method. 826 827 Input Parameter: 828 . snes - the SNES context 829 830 Application Interface Routine: SNESPrintHelp() 831 */ 832 #undef __FUNC__ 833 #define __FUNC__ "SNESPrintHelp_EQ_LS" 834 static int SNESPrintHelp_EQ_LS(SNES snes,char *p) 835 { 836 SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 837 int ierr; 838 839 PetscFunctionBegin; 840 ierr = (*PetscHelpPrintf)(snes->comm," method SNESEQLS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr); 841 ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr); 842 ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr); 843 ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr); 844 ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr); 845 PetscFunctionReturn(0); 846 } 847 /* -------------------------------------------------------------------------- */ 848 /* 849 SNESView_EQ_LS - Prints info from the SNESEQLS data structure. 850 851 Input Parameters: 852 . SNES - the SNES context 853 . viewer - visualization context 854 855 Application Interface Routine: SNESView() 856 */ 857 #undef __FUNC__ 858 #define __FUNC__ "SNESView_EQ_LS" 859 static int SNESView_EQ_LS(SNES snes,Viewer viewer) 860 { 861 SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 862 char *cstr; 863 int ierr; 864 PetscTruth isascii; 865 866 PetscFunctionBegin; 867 ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 868 if (isascii) { 869 if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 870 else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 871 else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 872 else cstr = "unknown"; 873 ierr = ViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 874 ierr = ViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 875 } else { 876 SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 877 } 878 PetscFunctionReturn(0); 879 } 880 /* -------------------------------------------------------------------------- */ 881 /* 882 SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method. 883 884 Input Parameter: 885 . snes - the SNES context 886 887 Application Interface Routine: SNESSetFromOptions() 888 */ 889 #undef __FUNC__ 890 #define __FUNC__ "SNESSetFromOptions_EQ_LS" 891 static int SNESSetFromOptions_EQ_LS(SNES snes) 892 { 893 SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data; 894 char ver[16]; 895 double tmp; 896 int ierr; 897 PetscTruth flg; 898 899 PetscFunctionBegin; 900 ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp,&flg);CHKERRQ(ierr); 901 if (flg) { 902 ls->alpha = tmp; 903 } 904 ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp,&flg);CHKERRQ(ierr); 905 if (flg) { 906 ls->maxstep = tmp; 907 } 908 ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp,&flg);CHKERRQ(ierr); 909 if (flg) { 910 ls->steptol = tmp; 911 } 912 ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16,&flg);CHKERRQ(ierr); 913 if (flg) { 914 PetscTruth isbasic,isnonorms,isquad,iscubic; 915 916 ierr = PetscStrcmp(ver,"basic",&isbasic);CHKERRQ(ierr); 917 ierr = PetscStrcmp(ver,"basicnonorms",&isnonorms);CHKERRQ(ierr); 918 ierr = PetscStrcmp(ver,"quadratic",&isquad);CHKERRQ(ierr); 919 ierr = PetscStrcmp(ver,"cubic",&iscubic);CHKERRQ(ierr); 920 921 if (isbasic) { 922 ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 923 } else if (isnonorms) { 924 ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 925 } else if (isquad) { 926 ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 927 } else if (iscubic) { 928 ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 929 } 930 else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");} 931 } 932 PetscFunctionReturn(0); 933 } 934 /* -------------------------------------------------------------------------- */ 935 /* 936 SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method, 937 SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver 938 context, SNES, that was created within SNESCreate(). 939 940 941 Input Parameter: 942 . snes - the SNES context 943 944 Application Interface Routine: SNESCreate() 945 */ 946 EXTERN_C_BEGIN 947 #undef __FUNC__ 948 #define __FUNC__ "SNESCreate_EQ_LS" 949 int SNESCreate_EQ_LS(SNES snes) 950 { 951 int ierr; 952 SNES_EQ_LS *neP; 953 954 PetscFunctionBegin; 955 if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 956 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 957 } 958 959 snes->setup = SNESSetUp_EQ_LS; 960 snes->solve = SNESSolve_EQ_LS; 961 snes->destroy = SNESDestroy_EQ_LS; 962 snes->converged = SNESConverged_EQ_LS; 963 snes->printhelp = SNESPrintHelp_EQ_LS; 964 snes->setfromoptions = SNESSetFromOptions_EQ_LS; 965 snes->view = SNESView_EQ_LS; 966 snes->nwork = 0; 967 968 neP = PetscNew(SNES_EQ_LS);CHKPTRQ(neP); 969 PLogObjectMemory(snes,sizeof(SNES_EQ_LS)); 970 snes->data = (void*)neP; 971 neP->alpha = 1.e-4; 972 neP->maxstep = 1.e8; 973 neP->steptol = 1.e-12; 974 neP->LineSearch = SNESCubicLineSearch; 975 neP->lsP = PETSC_NULL; 976 neP->CheckStep = PETSC_NULL; 977 neP->checkP = PETSC_NULL; 978 979 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS", 980 (void*)SNESSetLineSearch_LS);CHKERRQ(ierr); 981 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS", 982 (void*)SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 983 984 PetscFunctionReturn(0); 985 } 986 EXTERN_C_END 987 988 989 990