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