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