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