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