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