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