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