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) { 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 = W; snes->vec_sol_always = X; W = 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 316 . w - work vector 317 - fnorm - 2-norm of f 318 319 Output Parameters: 320 + g - residual evaluated at new iterate y 321 . w - new iterate 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 changed_w = PETSC_FALSE,changed_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 = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr); /* w <- x - y */ 348 if (neP->postcheckstep) { 349 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 350 } 351 if (changed_y) { 352 ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr); /* w <- x - y */ 353 } 354 ierr = SNESComputeFunction(snes,w,g); 355 if (PetscExceptionValue(ierr)) { 356 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 357 } 358 CHKERRQ(ierr); 359 360 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 361 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 362 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 363 PetscFunctionReturn(0); 364 } 365 /* -------------------------------------------------------------------------- */ 366 #undef __FUNCT__ 367 #define __FUNCT__ "SNESNoLineSearchNoNorms" 368 369 /*@C 370 SNESNoLineSearchNoNorms - This routine is not a line search at 371 all; it simply uses the full Newton step. This version does not 372 even compute the norm of the function or search direction; this 373 is intended only when you know the full step is fine and are 374 not checking for convergence of the nonlinear iteration (for 375 example, you are running always for a fixed number of Newton steps). 376 377 Collective on SNES and Vec 378 379 Input Parameters: 380 + snes - nonlinear context 381 . lsctx - optional context for line search (not used here) 382 . x - current iterate 383 . f - residual evaluated at x 384 . y - search direction 385 . w - work vector 386 - fnorm - 2-norm of f 387 388 Output Parameters: 389 + g - residual evaluated at new iterate y 390 . w - new iterate 391 . gnorm - not changed 392 . ynorm - not changed 393 - flag - set to PETSC_TRUE indicating a successful line search 394 395 Options Database Key: 396 . -snes_ls basicnonorms - Activates SNESNoLineSearchNoNorms() 397 398 Notes: 399 SNESNoLineSearchNoNorms() must be used in conjunction with 400 either the options 401 $ -snes_no_convergence_test -snes_max_it <its> 402 or alternatively a user-defined custom test set via 403 SNESSetConvergenceTest(); or a -snes_max_it of 1, 404 otherwise, the SNES solver will generate an error. 405 406 During the final iteration this will not evaluate the function at 407 the solution point. This is to save a function evaluation while 408 using pseudo-timestepping. 409 410 The residual norms printed by monitoring routines such as 411 SNESDefaultMonitor() (as activated via -snes_monitor) will not be 412 correct, since they are not computed. 413 414 Level: advanced 415 416 .keywords: SNES, nonlinear, line search, cubic 417 418 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), 419 SNESSetLineSearch(), SNESNoLineSearch() 420 @*/ 421 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) 422 { 423 PetscErrorCode ierr; 424 PetscScalar mone = -1.0; 425 SNES_LS *neP = (SNES_LS*)snes->data; 426 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 427 428 PetscFunctionBegin; 429 *flag = PETSC_TRUE; 430 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 431 ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr); /* w <- x - y */ 432 if (neP->postcheckstep) { 433 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 434 } 435 if (changed_y) { 436 ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr); /* w <- x - y */ 437 } 438 439 /* don't evaluate function the last time through */ 440 if (snes->iter < snes->max_its-1) { 441 ierr = SNESComputeFunction(snes,w,g); 442 if (PetscExceptionValue(ierr)) { 443 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 444 } 445 CHKERRQ(ierr); 446 } 447 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 448 PetscFunctionReturn(0); 449 } 450 /* -------------------------------------------------------------------------- */ 451 #undef __FUNCT__ 452 #define __FUNCT__ "SNESCubicLineSearch" 453 /*@C 454 SNESCubicLineSearch - Performs a cubic line search (default line search method). 455 456 Collective on SNES 457 458 Input Parameters: 459 + snes - nonlinear context 460 . lsctx - optional context for line search (not used here) 461 . x - current iterate 462 . f - residual evaluated at x 463 . y - search direction 464 . w - work vector 465 - fnorm - 2-norm of f 466 467 Output Parameters: 468 + g - residual evaluated at new iterate y 469 . w - new iterate 470 . gnorm - 2-norm of g 471 . ynorm - 2-norm of search length 472 - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 473 474 Options Database Key: 475 $ -snes_ls cubic - Activates SNESCubicLineSearch() 476 477 Notes: 478 This line search is taken from "Numerical Methods for Unconstrained 479 Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 480 481 Level: advanced 482 483 .keywords: SNES, nonlinear, line search, cubic 484 485 .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 486 @*/ 487 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) 488 { 489 /* 490 Note that for line search purposes we work with with the related 491 minimization problem: 492 min z(x): R^n -> R, 493 where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 494 */ 495 496 PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 497 PetscReal maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg; 498 #if defined(PETSC_USE_COMPLEX) 499 PetscScalar cinitslope,clambda; 500 #endif 501 PetscErrorCode ierr; 502 PetscInt count; 503 SNES_LS *neP = (SNES_LS*)snes->data; 504 PetscScalar mone = -1.0,scale; 505 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 506 507 PetscFunctionBegin; 508 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 509 *flag = PETSC_TRUE; 510 alpha = neP->alpha; 511 maxstep = neP->maxstep; 512 steptol = neP->steptol; 513 514 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 515 if (!*ynorm) { 516 ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Search direction and size is 0\n"));CHKERRQ(ierr); 517 *gnorm = fnorm; 518 ierr = VecCopy(x,w);CHKERRQ(ierr); 519 ierr = VecCopy(f,g);CHKERRQ(ierr); 520 *flag = PETSC_FALSE; 521 goto theend1; 522 } 523 if (*ynorm > maxstep) { /* Step too big, so scale back */ 524 scale = maxstep/(*ynorm); 525 #if defined(PETSC_USE_COMPLEX) 526 ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm));CHKERRQ(ierr); 527 #else 528 ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm));CHKERRQ(ierr); 529 #endif 530 ierr = VecScale(&scale,y);CHKERRQ(ierr); 531 *ynorm = maxstep; 532 } 533 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 534 minlambda = steptol/rellength; 535 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 536 #if defined(PETSC_USE_COMPLEX) 537 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 538 initslope = PetscRealPart(cinitslope); 539 #else 540 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 541 #endif 542 if (initslope > 0.0) initslope = -initslope; 543 if (initslope == 0.0) initslope = -1.0; 544 545 ierr = VecCopy(y,w);CHKERRQ(ierr); 546 ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 547 if (snes->nfuncs >= snes->max_funcs) { 548 ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr); 549 *flag = PETSC_FALSE; 550 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 551 goto theend1; 552 } 553 ierr = SNESComputeFunction(snes,w,g); 554 if (PetscExceptionValue(ierr)) { 555 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 556 } 557 CHKERRQ(ierr); 558 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 559 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 560 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 561 ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Using full step\n"));CHKERRQ(ierr); 562 goto theend1; 563 } 564 565 /* Fit points with quadratic */ 566 lambda = 1.0; 567 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 568 lambdaprev = lambda; 569 gnormprev = *gnorm; 570 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 571 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 572 else lambda = lambdatemp; 573 ierr = VecCopy(x,w);CHKERRQ(ierr); 574 lambdaneg = -lambda; 575 #if defined(PETSC_USE_COMPLEX) 576 clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 577 #else 578 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 579 #endif 580 if (snes->nfuncs >= snes->max_funcs) { 581 ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n"));CHKERRQ(ierr); 582 *flag = PETSC_FALSE; 583 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 584 goto theend1; 585 } 586 ierr = SNESComputeFunction(snes,w,g); 587 if (PetscExceptionValue(ierr)) { 588 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 589 } 590 CHKERRQ(ierr); 591 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 592 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 593 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 594 ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr); 595 goto theend1; 596 } 597 598 /* Fit points with cubic */ 599 count = 1; 600 while (count < 10000) { 601 if (lambda <= minlambda) { /* bad luck; use full step */ 602 ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Unable to find good step length! %D \n",count));CHKERRQ(ierr); 603 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); 604 *flag = PETSC_FALSE; 605 break; 606 } 607 t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 608 t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 609 a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 610 b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 611 d = b*b - 3*a*initslope; 612 if (d < 0.0) d = 0.0; 613 if (a == 0.0) { 614 lambdatemp = -initslope/(2.0*b); 615 } else { 616 lambdatemp = (-b + sqrt(d))/(3.0*a); 617 } 618 lambdaprev = lambda; 619 gnormprev = *gnorm; 620 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 621 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 622 else lambda = lambdatemp; 623 ierr = VecCopy(x,w);CHKERRQ(ierr); 624 lambdaneg = -lambda; 625 #if defined(PETSC_USE_COMPLEX) 626 clambda = lambdaneg; 627 ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 628 #else 629 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 630 #endif 631 if (snes->nfuncs >= snes->max_funcs) { 632 ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr); 633 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); 634 *flag = PETSC_FALSE; 635 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 636 break; 637 } 638 ierr = SNESComputeFunction(snes,w,g); 639 if (PetscExceptionValue(ierr)) { 640 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 641 } 642 CHKERRQ(ierr); 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 = PetscLogInfo((snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr); 647 break; 648 } else { 649 ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda, lambda=%18.16e\n",lambda));CHKERRQ(ierr); 650 } 651 count++; 652 } 653 if (count >= 10000) { 654 SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation"); 655 } 656 theend1: 657 /* Optional user-defined check for line search step validity */ 658 if (neP->postcheckstep && *flag) { 659 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 660 if (changed_y) { 661 ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr); 662 } 663 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 664 ierr = SNESComputeFunction(snes,w,g); 665 if (PetscExceptionValue(ierr)) { 666 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 667 } 668 CHKERRQ(ierr); 669 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 670 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 671 ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr); 672 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 673 ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr); 674 } 675 } 676 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 677 PetscFunctionReturn(0); 678 } 679 /* -------------------------------------------------------------------------- */ 680 #undef __FUNCT__ 681 #define __FUNCT__ "SNESQuadraticLineSearch" 682 /*@C 683 SNESQuadraticLineSearch - Performs a quadratic line search. 684 685 Collective on SNES and Vec 686 687 Input Parameters: 688 + snes - the SNES context 689 . lsctx - optional context for line search (not used here) 690 . x - current iterate 691 . f - residual evaluated at x 692 . y - search direction 693 . w - work vector 694 - fnorm - 2-norm of f 695 696 Output Parameters: 697 + g - residual evaluated at new iterate y 698 . w - new iterate 699 . gnorm - 2-norm of g 700 . ynorm - 2-norm of search length 701 - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 702 703 Options Database Key: 704 . -snes_ls quadratic - Activates SNESQuadraticLineSearch() 705 706 Notes: 707 Use SNESSetLineSearch() to set this routine within the SNESLS method. 708 709 Level: advanced 710 711 .keywords: SNES, nonlinear, quadratic, line search 712 713 .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms() 714 @*/ 715 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) 716 { 717 /* 718 Note that for line search purposes we work with with the related 719 minimization problem: 720 min z(x): R^n -> R, 721 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 722 */ 723 PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength; 724 #if defined(PETSC_USE_COMPLEX) 725 PetscScalar cinitslope,clambda; 726 #endif 727 PetscErrorCode ierr; 728 PetscInt count; 729 SNES_LS *neP = (SNES_LS*)snes->data; 730 PetscScalar mone = -1.0,scale; 731 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 732 733 PetscFunctionBegin; 734 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 735 *flag = PETSC_TRUE; 736 alpha = neP->alpha; 737 maxstep = neP->maxstep; 738 steptol = neP->steptol; 739 740 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 741 if (*ynorm == 0.0) { 742 ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"));CHKERRQ(ierr); 743 *gnorm = fnorm; 744 ierr = VecCopy(x,y);CHKERRQ(ierr); 745 ierr = VecCopy(f,g);CHKERRQ(ierr); 746 *flag = PETSC_FALSE; 747 goto theend2; 748 } 749 if (*ynorm > maxstep) { /* Step too big, so scale back */ 750 scale = maxstep/(*ynorm); 751 ierr = VecScale(&scale,y);CHKERRQ(ierr); 752 *ynorm = maxstep; 753 } 754 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 755 minlambda = steptol/rellength; 756 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 757 #if defined(PETSC_USE_COMPLEX) 758 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 759 initslope = PetscRealPart(cinitslope); 760 #else 761 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 762 #endif 763 if (initslope > 0.0) initslope = -initslope; 764 if (initslope == 0.0) initslope = -1.0; 765 766 ierr = VecCopy(y,w);CHKERRQ(ierr); 767 ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr); 768 if (snes->nfuncs >= snes->max_funcs) { 769 ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr); 770 ierr = VecCopy(w,y);CHKERRQ(ierr); 771 *flag = PETSC_FALSE; 772 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 773 goto theend2; 774 } 775 ierr = SNESComputeFunction(snes,w,g); 776 if (PetscExceptionValue(ierr)) { 777 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 778 } 779 CHKERRQ(ierr); 780 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 781 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 782 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 783 ierr = VecCopy(w,y);CHKERRQ(ierr); 784 ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch: Using full step\n"));CHKERRQ(ierr); 785 goto theend2; 786 } 787 788 /* Fit points with quadratic */ 789 lambda = 1.0; 790 count = 1; 791 while (PETSC_TRUE) { 792 if (lambda <= minlambda) { /* bad luck; use full step */ 793 ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Unable to find good step length! %D \n",count));CHKERRQ(ierr); 794 ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr); 795 ierr = VecCopy(x,y);CHKERRQ(ierr); 796 *flag = PETSC_FALSE; 797 break; 798 } 799 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 800 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 801 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 802 else lambda = lambdatemp; 803 ierr = VecCopy(x,w);CHKERRQ(ierr); 804 lambdaneg = -lambda; 805 #if defined(PETSC_USE_COMPLEX) 806 clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr); 807 #else 808 ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr); 809 #endif 810 if (snes->nfuncs >= snes->max_funcs) { 811 ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr); 812 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); 813 ierr = VecCopy(w,y);CHKERRQ(ierr); 814 *flag = PETSC_FALSE; 815 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 816 break; 817 } 818 ierr = SNESComputeFunction(snes,w,g); 819 if (PetscExceptionValue(ierr)) { 820 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 821 } 822 CHKERRQ(ierr); 823 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 824 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 825 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 826 ierr = VecCopy(w,y);CHKERRQ(ierr); 827 ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda));CHKERRQ(ierr); 828 break; 829 } 830 count++; 831 } 832 theend2: 833 /* Optional user-defined check for line search step validity */ 834 if (neP->postcheckstep) { 835 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 836 if (changed_y) { 837 ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr); 838 } 839 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 840 ierr = SNESComputeFunction(snes,w,g); 841 if (PetscExceptionValue(ierr)) { 842 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 843 } 844 CHKERRQ(ierr); 845 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 846 ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr); 847 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 848 ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr); 849 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 850 } 851 } 852 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 853 PetscFunctionReturn(0); 854 } 855 856 /* -------------------------------------------------------------------------- */ 857 #undef __FUNCT__ 858 #define __FUNCT__ "SNESSetLineSearch" 859 /*@C 860 SNESSetLineSearch - Sets the line search routine to be used 861 by the method SNESLS. 862 863 Input Parameters: 864 + snes - nonlinear context obtained from SNESCreate() 865 . lsctx - optional user-defined context for use by line search 866 - func - pointer to int function 867 868 Collective on SNES 869 870 Available Routines: 871 + SNESCubicLineSearch() - default line search 872 . SNESQuadraticLineSearch() - quadratic line search 873 . SNESNoLineSearch() - the full Newton step (actually not a line search) 874 - SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 875 876 Options Database Keys: 877 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 878 . -snes_ls_alpha <alpha> - Sets alpha 879 . -snes_ls_maxstep <max> - Sets maxstep 880 - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 881 will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 882 the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 883 884 Calling sequence of func: 885 .vb 886 func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 887 PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 888 .ve 889 890 Input parameters for func: 891 + snes - nonlinear context 892 . lsctx - optional user-defined context for line search 893 . x - current iterate 894 . f - residual evaluated at x 895 . y - search direction 896 . w - work vector 897 - fnorm - 2-norm of f 898 899 Output parameters for func: 900 + g - residual evaluated at new iterate y 901 . w - new iterate 902 . gnorm - 2-norm of g 903 . ynorm - 2-norm of search length 904 - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 905 906 Level: advanced 907 908 .keywords: SNES, nonlinear, set, line search, routine 909 910 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), 911 SNESLineSearchSetPostCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams() 912 @*/ 913 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLineSearch(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx) 914 { 915 PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*); 916 917 PetscFunctionBegin; 918 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr); 919 if (f) { 920 ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 921 } 922 PetscFunctionReturn(0); 923 } 924 925 typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/ 926 /* -------------------------------------------------------------------------- */ 927 EXTERN_C_BEGIN 928 #undef __FUNCT__ 929 #define __FUNCT__ "SNESSetLineSearch_LS" 930 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLineSearch_LS(SNES snes,FCN2 func,void *lsctx) 931 { 932 PetscFunctionBegin; 933 ((SNES_LS *)(snes->data))->LineSearch = func; 934 ((SNES_LS *)(snes->data))->lsP = lsctx; 935 PetscFunctionReturn(0); 936 } 937 EXTERN_C_END 938 /* -------------------------------------------------------------------------- */ 939 #undef __FUNCT__ 940 #define __FUNCT__ "SNESLineSearchSetPostCheck" 941 /*@C 942 SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed 943 by the line search routine in the Newton-based method SNESLS. 944 945 Input Parameters: 946 + snes - nonlinear context obtained from SNESCreate() 947 . func - pointer to function 948 - checkctx - optional user-defined context for use by step checking routine 949 950 Collective on SNES 951 952 Calling sequence of func: 953 .vb 954 int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w) 955 .ve 956 where func returns an error code of 0 on success and a nonzero 957 on failure. 958 959 Input parameters for func: 960 + snes - nonlinear context 961 . checkctx - optional user-defined context for use by step checking routine 962 . x - previous iterate 963 . y - new search direction and length 964 - w - current candidate iterate 965 966 Output parameters for func: 967 + y - search direction (possibly changed) 968 . w - current iterate (possibly modified) 969 . changed_y - indicates search direction was changed by this routine 970 - changed_w - indicates current iterate was changed by this routine 971 972 Level: advanced 973 974 Notes: All line searches accept the new iterate computed by the line search checking routine. 975 976 Only one of changed_y and changed_w can be PETSC_TRUE 977 978 On input w = x + y 979 980 SNESNoLineSearch() and SNESNoLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control 981 to the checking routine, and then (3) compute the corresponding nonlinear 982 function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 983 984 SNESQuadraticLineSearch() and SNESCubicLineSearch() (1) compute a candidate iterate u_{i+1} as well as a 985 candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 986 routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 987 were made to the candidate iterate in the checking routine (as indicated 988 by flag=PETSC_TRUE). The overhead of this extra function re-evaluation can be 989 very costly, so use this feature with caution! 990 991 .keywords: SNES, nonlinear, set, line search check, step check, routine 992 993 .seealso: SNESSetLineSearch() 994 @*/ 995 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx) 996 { 997 PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*); 998 999 PetscFunctionBegin; 1000 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 1001 if (f) { 1002 ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 1003 } 1004 PetscFunctionReturn(0); 1005 } 1006 /* -------------------------------------------------------------------------- */ 1007 typedef PetscErrorCode (*FCN)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/ 1008 EXTERN_C_BEGIN 1009 #undef __FUNCT__ 1010 #define __FUNCT__ "SNESLineSearchSetPostCheck_LS" 1011 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN func,void *checkctx) 1012 { 1013 PetscFunctionBegin; 1014 ((SNES_LS *)(snes->data))->postcheckstep = func; 1015 ((SNES_LS *)(snes->data))->postcheck = checkctx; 1016 PetscFunctionReturn(0); 1017 } 1018 EXTERN_C_END 1019 /* -------------------------------------------------------------------------- */ 1020 /* 1021 SNESPrintHelp_LS - Prints all options for the SNES_LS method. 1022 1023 Input Parameter: 1024 . snes - the SNES context 1025 1026 Application Interface Routine: SNESPrintHelp() 1027 */ 1028 #undef __FUNCT__ 1029 #define __FUNCT__ "SNESPrintHelp_LS" 1030 static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p) 1031 { 1032 SNES_LS *ls = (SNES_LS *)snes->data; 1033 1034 PetscFunctionBegin; 1035 (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n"); 1036 (*PetscHelpPrintf)(snes->comm," %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p); 1037 (*PetscHelpPrintf)(snes->comm," %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha); 1038 (*PetscHelpPrintf)(snes->comm," %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep); 1039 (*PetscHelpPrintf)(snes->comm," %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol); 1040 PetscFunctionReturn(0); 1041 } 1042 1043 /* 1044 SNESView_LS - Prints info from the SNESLS data structure. 1045 1046 Input Parameters: 1047 . SNES - the SNES context 1048 . viewer - visualization context 1049 1050 Application Interface Routine: SNESView() 1051 */ 1052 #undef __FUNCT__ 1053 #define __FUNCT__ "SNESView_LS" 1054 static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 1055 { 1056 SNES_LS *ls = (SNES_LS *)snes->data; 1057 const char *cstr; 1058 PetscErrorCode ierr; 1059 PetscTruth iascii; 1060 1061 PetscFunctionBegin; 1062 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1063 if (iascii) { 1064 if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch"; 1065 else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch"; 1066 else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch"; 1067 else cstr = "unknown"; 1068 ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1069 ierr = PetscViewerASCIIPrintf(viewer," alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 1070 } else { 1071 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 1072 } 1073 PetscFunctionReturn(0); 1074 } 1075 /* -------------------------------------------------------------------------- */ 1076 /* 1077 SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 1078 1079 Input Parameter: 1080 . snes - the SNES context 1081 1082 Application Interface Routine: SNESSetFromOptions() 1083 */ 1084 #undef __FUNCT__ 1085 #define __FUNCT__ "SNESSetFromOptions_LS" 1086 static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 1087 { 1088 SNES_LS *ls = (SNES_LS *)snes->data; 1089 const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1090 PetscErrorCode ierr; 1091 PetscInt indx; 1092 PetscTruth flg; 1093 1094 PetscFunctionBegin; 1095 ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 1096 ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 1097 ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 1098 ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr); 1099 1100 ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 1101 if (flg) { 1102 switch (indx) { 1103 case 0: 1104 ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr); 1105 break; 1106 case 1: 1107 ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1108 break; 1109 case 2: 1110 ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr); 1111 break; 1112 case 3: 1113 ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr); 1114 break; 1115 } 1116 } 1117 ierr = PetscOptionsTail();CHKERRQ(ierr); 1118 PetscFunctionReturn(0); 1119 } 1120 /* -------------------------------------------------------------------------- */ 1121 /*MC 1122 SNESLS - Newton based nonlinear solver that uses a line search 1123 1124 Options Database: 1125 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1126 . -snes_ls_alpha <alpha> - Sets alpha 1127 . -snes_ls_maxstep <max> - Sets maxstep 1128 - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 1129 will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 1130 the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 1131 1132 Notes: This is the default nonlinear solver in SNES 1133 1134 Level: beginner 1135 1136 .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESSetLineSearch(), 1137 SNESLineSearchSetPostCheck(), SNESNoLineSearch(), SNESCubicLineSearch(), SNESQuadraticLineSearch(), 1138 SNESSetLineSearch(), SNESNoLineSearchNoNorms() 1139 1140 M*/ 1141 EXTERN_C_BEGIN 1142 #undef __FUNCT__ 1143 #define __FUNCT__ "SNESCreate_LS" 1144 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes) 1145 { 1146 PetscErrorCode ierr; 1147 SNES_LS *neP; 1148 1149 PetscFunctionBegin; 1150 snes->setup = SNESSetUp_LS; 1151 snes->solve = SNESSolve_LS; 1152 snes->destroy = SNESDestroy_LS; 1153 snes->converged = SNESConverged_LS; 1154 snes->printhelp = SNESPrintHelp_LS; 1155 snes->setfromoptions = SNESSetFromOptions_LS; 1156 snes->view = SNESView_LS; 1157 snes->nwork = 0; 1158 1159 ierr = PetscNew(SNES_LS,&neP);CHKERRQ(ierr); 1160 ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr); 1161 snes->data = (void*)neP; 1162 neP->alpha = 1.e-4; 1163 neP->maxstep = 1.e8; 1164 neP->steptol = 1.e-12; 1165 neP->LineSearch = SNESCubicLineSearch; 1166 neP->lsP = PETSC_NULL; 1167 neP->postcheckstep = PETSC_NULL; 1168 neP->postcheck = PETSC_NULL; 1169 neP->precheckstep = PETSC_NULL; 1170 neP->precheck = PETSC_NULL; 1171 1172 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr); 1173 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr); 1174 1175 PetscFunctionReturn(0); 1176 } 1177 EXTERN_C_END 1178 1179 1180 1181