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