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 int SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin) 13 { 14 PetscReal a1; 15 int 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 int SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2) 51 { 52 PetscReal a1,a2; 53 int 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 int 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 int SNESSetUp_LS(SNES snes) 261 { 262 int 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 int SNESDestroy_LS(SNES snes) 284 { 285 int 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 int 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 int 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 int 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 int 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 int 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 int ierr,count; 480 SNES_LS *neP = (SNES_LS*)snes->data; 481 PetscScalar mone = -1.0,scale; 482 PetscTruth change_y = PETSC_FALSE; 483 484 PetscFunctionBegin; 485 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 486 *flag = 0; 487 alpha = neP->alpha; 488 maxstep = neP->maxstep; 489 steptol = neP->steptol; 490 491 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 492 if (*ynorm == 0.0) { 493 PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size is 0\n"); 494 *gnorm = fnorm; 495 ierr = VecCopy(x,y);CHKERRQ(ierr); 496 ierr = VecCopy(f,g);CHKERRQ(ierr); 497 *flag = -1; 498 goto theend1; 499 } 500 if (*ynorm > maxstep) { /* Step too big, so scale back */ 501 scale = maxstep/(*ynorm); 502 #if defined(PETSC_USE_COMPLEX) 503 PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm); 504 #else 505 PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm); 506 #endif 507 ierr = VecScale(&scale,y);CHKERRQ(ierr); 508 *ynorm = maxstep; 509 } 510 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 511 minlambda = steptol/rellength; 512 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 513 #if defined(PETSC_USE_COMPLEX) 514 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 515 initslope = PetscRealPart(cinitslope); 516 #else 517 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 518 #endif 519 if (initslope > 0.0) initslope = -initslope; 520 if (initslope == 0.0) initslope = -1.0; 521 522 ierr = VecCopy(y,w);CHKERRQ(ierr); 523 ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 524 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 525 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 526 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 527 ierr = VecCopy(w,y);CHKERRQ(ierr); 528 PetscLogInfo(snes,"SNESCubicLineSearch: Using full step\n"); 529 goto theend1; 530 } 531 532 /* Fit points with quadratic */ 533 lambda = 1.0; 534 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 535 lambdaprev = lambda; 536 gnormprev = *gnorm; 537 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 538 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 539 else lambda = lambdatemp; 540 ierr = VecCopy(x,w);CHKERRQ(ierr); 541 lambdaneg = -lambda; 542 #if defined(PETSC_USE_COMPLEX) 543 clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 544 #else 545 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 546 #endif 547 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 548 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 549 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 550 ierr = VecCopy(w,y);CHKERRQ(ierr); 551 PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda); 552 goto theend1; 553 } 554 555 /* Fit points with cubic */ 556 count = 1; 557 while (count < 10000) { 558 if (lambda <= minlambda) { /* bad luck; use full step */ 559 PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count); 560 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); 561 ierr = VecCopy(x,y);CHKERRQ(ierr); 562 *flag = -1; break; 563 } 564 t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 565 t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 566 a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 567 b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 568 d = b*b - 3*a*initslope; 569 if (d < 0.0) d = 0.0; 570 if (a == 0.0) { 571 lambdatemp = -initslope/(2.0*b); 572 } else { 573 lambdatemp = (-b + sqrt(d))/(3.0*a); 574 } 575 lambdaprev = lambda; 576 gnormprev = *gnorm; 577 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 578 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 579 else lambda = lambdatemp; 580 ierr = VecCopy(x,w);CHKERRQ(ierr); 581 lambdaneg = -lambda; 582 #if defined(PETSC_USE_COMPLEX) 583 clambda = lambdaneg; 584 ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 585 #else 586 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 587 #endif 588 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 589 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 590 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 591 ierr = VecCopy(w,y);CHKERRQ(ierr); 592 PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda); 593 goto theend1; 594 } else { 595 PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%18.16e\n",lambda); 596 } 597 count++; 598 } 599 if (count >= 10000) { 600 SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation"); 601 } 602 theend1: 603 /* Optional user-defined check for line search step validity */ 604 if (neP->CheckStep) { 605 ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 606 if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 607 ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 608 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 609 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 610 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 611 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 612 } 613 } 614 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 615 PetscFunctionReturn(0); 616 } 617 /* -------------------------------------------------------------------------- */ 618 #undef __FUNCT__ 619 #define __FUNCT__ "SNESQuadraticLineSearch" 620 /*@C 621 SNESQuadraticLineSearch - Performs a quadratic line search. 622 623 Collective on SNES and Vec 624 625 Input Parameters: 626 + snes - the SNES context 627 . lsctx - optional context for line search (not used here) 628 . x - current iterate 629 . f - residual evaluated at x 630 . y - search direction (contains new iterate on output) 631 . w - work vector 632 - fnorm - 2-norm of f 633 634 Output Parameters: 635 + g - residual evaluated at new iterate y 636 . y - new iterate (contains search direction on input) 637 . gnorm - 2-norm of g 638 . ynorm - 2-norm of search length 639 - flag - 0 if line search succeeds; -1 on failure. 640 641 Options Database Key: 642 . -snes_ls quadratic - Activates SNESQuadraticLineSearch() 643 644 Notes: 645 Use SNESSetLineSearch() to set this routine within the SNESLS method. 646 647 Level: advanced 648 649 .keywords: SNES, nonlinear, quadratic, line search 650 651 .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 652 @*/ 653 int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag) 654 { 655 /* 656 Note that for line search purposes we work with with the related 657 minimization problem: 658 min z(x): R^n -> R, 659 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 660 */ 661 PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength; 662 #if defined(PETSC_USE_COMPLEX) 663 PetscScalar cinitslope,clambda; 664 #endif 665 int ierr,count; 666 SNES_LS *neP = (SNES_LS*)snes->data; 667 PetscScalar mone = -1.0,scale; 668 PetscTruth change_y = PETSC_FALSE; 669 670 PetscFunctionBegin; 671 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 672 *flag = 0; 673 alpha = neP->alpha; 674 maxstep = neP->maxstep; 675 steptol = neP->steptol; 676 677 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 678 if (*ynorm == 0.0) { 679 PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"); 680 *gnorm = fnorm; 681 ierr = VecCopy(x,y);CHKERRQ(ierr); 682 ierr = VecCopy(f,g);CHKERRQ(ierr); 683 *flag = -1; 684 goto theend2; 685 } 686 if (*ynorm > maxstep) { /* Step too big, so scale back */ 687 scale = maxstep/(*ynorm); 688 ierr = VecScale(&scale,y);CHKERRQ(ierr); 689 *ynorm = maxstep; 690 } 691 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 692 minlambda = steptol/rellength; 693 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 694 #if defined(PETSC_USE_COMPLEX) 695 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 696 initslope = PetscRealPart(cinitslope); 697 #else 698 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 699 #endif 700 if (initslope > 0.0) initslope = -initslope; 701 if (initslope == 0.0) initslope = -1.0; 702 703 ierr = VecCopy(y,w);CHKERRQ(ierr); 704 ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 705 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 706 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 707 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 708 ierr = VecCopy(w,y);CHKERRQ(ierr); 709 PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n"); 710 goto theend2; 711 } 712 713 /* Fit points with quadratic */ 714 lambda = 1.0; 715 count = 1; 716 while (PETSC_TRUE) { 717 if (lambda <= minlambda) { /* bad luck; use full step */ 718 PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count); 719 PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope); 720 ierr = VecCopy(x,y);CHKERRQ(ierr); 721 *flag = -1; break; 722 } 723 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 724 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 725 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 726 else lambda = lambdatemp; 727 ierr = VecCopy(x,w);CHKERRQ(ierr); 728 lambdaneg = -lambda; 729 #if defined(PETSC_USE_COMPLEX) 730 clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 731 #else 732 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 733 #endif 734 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 735 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 736 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 737 ierr = VecCopy(w,y);CHKERRQ(ierr); 738 PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda); 739 break; 740 } 741 count++; 742 } 743 theend2: 744 /* Optional user-defined check for line search step validity */ 745 if (neP->CheckStep) { 746 ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr); 747 if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */ 748 ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); 749 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 750 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 751 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 752 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 753 } 754 } 755 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 756 PetscFunctionReturn(0); 757 } 758 759 /* -------------------------------------------------------------------------- */ 760 #undef __FUNCT__ 761 #define __FUNCT__ "SNESSetLineSearch" 762 /*@C 763 SNESSetLineSearch - Sets the line search routine to be used 764 by the method SNESLS. 765 766 Input Parameters: 767 + snes - nonlinear context obtained from SNESCreate() 768 . lsctx - optional user-defined context for use by line search 769 - func - pointer to int function 770 771 Collective on SNES 772 773 Available Routines: 774 + SNESCubicLineSearch() - default line search 775 . SNESQuadraticLineSearch() - quadratic line search 776 . SNESNoLineSearch() - the full Newton step (actually not a line search) 777 - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 778 779 Options Database Keys: 780 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 781 . -snes_ls_alpha <alpha> - Sets alpha 782 . -snes_ls_maxstep <max> - Sets maxstep 783 - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 784 will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 785 the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 786 787 Calling sequence of func: 788 .vb 789 func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 790 PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag) 791 .ve 792 793 Input parameters for func: 794 + snes - nonlinear context 795 . lsctx - optional user-defined context for line search 796 . x - current iterate 797 . f - residual evaluated at x 798 . y - search direction (contains new iterate on output) 799 . w - work vector 800 - fnorm - 2-norm of f 801 802 Output parameters for func: 803 + g - residual evaluated at new iterate y 804 . y - new iterate (contains search direction on input) 805 . gnorm - 2-norm of g 806 . ynorm - 2-norm of search length 807 - flag - set to 0 if the line search succeeds; a nonzero integer 808 on failure. 809 810 Level: advanced 811 812 .keywords: SNES, nonlinear, set, line search, routine 813 814 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), 815 SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 816 @*/ 817 int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx) 818 { 819 int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*); 820 821 PetscFunctionBegin; 822 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr); 823 if (f) { 824 ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 825 } 826 PetscFunctionReturn(0); 827 } 828 829 typedef int (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*); /* force argument to next function to not be extern C*/ 830 /* -------------------------------------------------------------------------- */ 831 EXTERN_C_BEGIN 832 #undef __FUNCT__ 833 #define __FUNCT__ "SNESSetLineSearch_LS" 834 int SNESSetLineSearch_LS(SNES snes,FCN2 func,void *lsctx) 835 { 836 PetscFunctionBegin; 837 ((SNES_LS *)(snes->data))->LineSearch = func; 838 ((SNES_LS *)(snes->data))->lsP = lsctx; 839 PetscFunctionReturn(0); 840 } 841 EXTERN_C_END 842 /* -------------------------------------------------------------------------- */ 843 #undef __FUNCT__ 844 #define __FUNCT__ "SNESSetLineSearchCheck" 845 /*@C 846 SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed 847 by the line search routine in the Newton-based method SNESLS. 848 849 Input Parameters: 850 + snes - nonlinear context obtained from SNESCreate() 851 . func - pointer to int function 852 - checkctx - optional user-defined context for use by step checking routine 853 854 Collective on SNES 855 856 Calling sequence of func: 857 .vb 858 int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag) 859 .ve 860 where func returns an error code of 0 on success and a nonzero 861 on failure. 862 863 Input parameters for func: 864 + snes - nonlinear context 865 . checkctx - optional user-defined context for use by step checking routine 866 - x - current candidate iterate 867 868 Output parameters for func: 869 + x - current iterate (possibly modified) 870 - flag - flag indicating whether x has been modified (either 871 PETSC_TRUE of PETSC_FALSE) 872 873 Level: advanced 874 875 Notes: 876 SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new 877 iterate computed by the line search checking routine. In particular, 878 these routines (1) compute a candidate iterate u_{i+1}, (2) pass control 879 to the checking routine, and then (3) compute the corresponding nonlinear 880 function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 881 882 SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the 883 new iterate computed by the line search checking routine. In particular, 884 these routines (1) compute a candidate iterate u_{i+1} as well as a 885 candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 886 routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 887 were made to the candidate iterate in the checking routine (as indicated 888 by flag=PETSC_TRUE). The overhead of this function re-evaluation can be 889 very costly, so use this feature with caution! 890 891 .keywords: SNES, nonlinear, set, line search check, step check, routine 892 893 .seealso: SNESSetLineSearch() 894 @*/ 895 int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx) 896 { 897 int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*); 898 899 PetscFunctionBegin; 900 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 901 if (f) { 902 ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 903 } 904 PetscFunctionReturn(0); 905 } 906 /* -------------------------------------------------------------------------- */ 907 typedef int (*FCN)(SNES,void*,Vec,PetscTruth*); /* force argument to next function to not be extern C*/ 908 EXTERN_C_BEGIN 909 #undef __FUNCT__ 910 #define __FUNCT__ "SNESSetLineSearchCheck_LS" 911 int SNESSetLineSearchCheck_LS(SNES snes,FCN func,void *checkctx) 912 { 913 PetscFunctionBegin; 914 ((SNES_LS *)(snes->data))->CheckStep = func; 915 ((SNES_LS *)(snes->data))->checkP = checkctx; 916 PetscFunctionReturn(0); 917 } 918 EXTERN_C_END 919 /* -------------------------------------------------------------------------- */ 920 /* 921 SNESPrintHelp_LS - Prints all options for the SNES_LS method. 922 923 Input Parameter: 924 . snes - the SNES context 925 926 Application Interface Routine: SNESPrintHelp() 927 */ 928 #undef __FUNCT__ 929 #define __FUNCT__ "SNESPrintHelp_LS" 930 static int SNESPrintHelp_LS(SNES snes,char *p) 931 { 932 SNES_LS *ls = (SNES_LS *)snes->data; 933 934 PetscFunctionBegin; 935 (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n"); 936 (*PetscHelpPrintf)(snes->comm," %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p); 937 (*PetscHelpPrintf)(snes->comm," %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha); 938 (*PetscHelpPrintf)(snes->comm," %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep); 939 (*PetscHelpPrintf)(snes->comm," %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol); 940 PetscFunctionReturn(0); 941 } 942 943 /* 944 SNESView_LS - Prints info from the SNESLS data structure. 945 946 Input Parameters: 947 . SNES - the SNES context 948 . viewer - visualization context 949 950 Application Interface Routine: SNESView() 951 */ 952 #undef __FUNCT__ 953 #define __FUNCT__ "SNESView_LS" 954 static int SNESView_LS(SNES snes,PetscViewer viewer) 955 { 956 SNES_LS *ls = (SNES_LS *)snes->data; 957 const char *cstr; 958 int ierr; 959 PetscTruth iascii; 960 961 PetscFunctionBegin; 962 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 963 if (iascii) { 964 if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 965 else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 966 else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 967 else cstr = "unknown"; 968 ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 969 ierr = PetscViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 970 } else { 971 SETERRQ1(1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 972 } 973 PetscFunctionReturn(0); 974 } 975 /* -------------------------------------------------------------------------- */ 976 /* 977 SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 978 979 Input Parameter: 980 . snes - the SNES context 981 982 Application Interface Routine: SNESSetFromOptions() 983 */ 984 #undef __FUNCT__ 985 #define __FUNCT__ "SNESSetFromOptions_LS" 986 static int SNESSetFromOptions_LS(SNES snes) 987 { 988 SNES_LS *ls = (SNES_LS *)snes->data; 989 const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 990 int ierr,indx; 991 PetscTruth flg; 992 993 PetscFunctionBegin; 994 ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 995 ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 996 ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 997 ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr); 998 999 ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 1000 if (flg) { 1001 switch (indx) { 1002 case 0: 1003 ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 1004 break; 1005 case 1: 1006 ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1007 break; 1008 case 2: 1009 ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 1010 break; 1011 case 3: 1012 ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 1013 break; 1014 } 1015 } 1016 ierr = PetscOptionsTail();CHKERRQ(ierr); 1017 PetscFunctionReturn(0); 1018 } 1019 /* -------------------------------------------------------------------------- */ 1020 /*MC 1021 SNESLS - Newton based nonlinear solver that uses a line search 1022 1023 Options Database: 1024 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1025 . -snes_ls_alpha <alpha> - Sets alpha 1026 . -snes_ls_maxstep <max> - Sets maxstep 1027 - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 1028 will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 1029 the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 1030 1031 Notes: This is the default nonlinear solver in SNES 1032 1033 Level: beginner 1034 1035 .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESSetLineSearch(), 1036 SNESSetLineSearchCheck(), SNESNoLineSearch(), SNESCubicLineSearch(), SNESQuadraticLineSearch(), 1037 SNESSetLineSearch(), SNESNoLineSearchNoNorms() 1038 1039 M*/ 1040 EXTERN_C_BEGIN 1041 #undef __FUNCT__ 1042 #define __FUNCT__ "SNESCreate_LS" 1043 int SNESCreate_LS(SNES snes) 1044 { 1045 int ierr; 1046 SNES_LS *neP; 1047 1048 PetscFunctionBegin; 1049 snes->setup = SNESSetUp_LS; 1050 snes->solve = SNESSolve_LS; 1051 snes->destroy = SNESDestroy_LS; 1052 snes->converged = SNESConverged_LS; 1053 snes->printhelp = SNESPrintHelp_LS; 1054 snes->setfromoptions = SNESSetFromOptions_LS; 1055 snes->view = SNESView_LS; 1056 snes->nwork = 0; 1057 1058 ierr = PetscNew(SNES_LS,&neP);CHKERRQ(ierr); 1059 PetscLogObjectMemory(snes,sizeof(SNES_LS)); 1060 snes->data = (void*)neP; 1061 neP->alpha = 1.e-4; 1062 neP->maxstep = 1.e8; 1063 neP->steptol = 1.e-12; 1064 neP->LineSearch = SNESCubicLineSearch; 1065 neP->lsP = PETSC_NULL; 1066 neP->CheckStep = PETSC_NULL; 1067 neP->checkP = PETSC_NULL; 1068 1069 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr); 1070 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr); 1071 1072 PetscFunctionReturn(0); 1073 } 1074 EXTERN_C_END 1075 1076 1077 1078