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