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