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 = PetscInfo1(0,"|| 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(W,PETSC_NULL);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 = PetscInfo1(0,"(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 57 PetscFunctionBegin; 58 ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 59 if (hastranspose) { 60 ierr = MatMult(A,X,W1);CHKERRQ(ierr); 61 ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 62 63 /* Compute || J^T W|| */ 64 ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 65 ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 66 ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 67 if (a1) { 68 ierr = PetscInfo1(0,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr); 69 } 70 } 71 PetscFunctionReturn(0); 72 } 73 74 /* -------------------------------------------------------------------- 75 76 This file implements a truncated Newton method with a line search, 77 for solving a system of nonlinear equations, using the KSP, Vec, 78 and Mat interfaces for linear solvers, vectors, and matrices, 79 respectively. 80 81 The following basic routines are required for each nonlinear solver: 82 SNESCreate_XXX() - Creates a nonlinear solver context 83 SNESSetFromOptions_XXX() - Sets runtime options 84 SNESSolve_XXX() - Solves the nonlinear system 85 SNESDestroy_XXX() - Destroys the nonlinear solver context 86 The suffix "_XXX" denotes a particular implementation, in this case 87 we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving 88 systems of nonlinear equations with a line search (LS) method. 89 These routines are actually called via the common user interface 90 routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 91 SNESDestroy(), so the application code interface remains identical 92 for all nonlinear solvers. 93 94 Another key routine is: 95 SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 96 by setting data structures and options. The interface routine SNESSetUp() 97 is not usually called directly by the user, but instead is called by 98 SNESSolve() if necessary. 99 100 Additional basic routines are: 101 SNESView_XXX() - Prints details of runtime options that 102 have actually been used. 103 These are called by application codes via the interface routines 104 SNESView(). 105 106 The various types of solvers (preconditioners, Krylov subspace methods, 107 nonlinear solvers, timesteppers) are all organized similarly, so the 108 above description applies to these categories also. 109 110 -------------------------------------------------------------------- */ 111 /* 112 SNESSolve_LS - Solves a nonlinear system with a truncated Newton 113 method with a line search. 114 115 Input Parameters: 116 . snes - the SNES context 117 118 Output Parameter: 119 . outits - number of iterations until termination 120 121 Application Interface Routine: SNESSolve() 122 123 Notes: 124 This implements essentially a truncated Newton method with a 125 line search. By default a cubic backtracking line search 126 is employed, as described in the text "Numerical Methods for 127 Unconstrained Optimization and Nonlinear Equations" by Dennis 128 and Schnabel. 129 */ 130 #undef __FUNCT__ 131 #define __FUNCT__ "SNESSolve_LS" 132 PetscErrorCode SNESSolve_LS(SNES snes) 133 { 134 SNES_LS *neP = (SNES_LS*)snes->data; 135 PetscErrorCode ierr; 136 PetscInt maxits,i,lits; 137 PetscTruth lssucceed; 138 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 139 PetscReal fnorm,gnorm,xnorm,ynorm; 140 Vec Y,X,F,G,W,TMP; 141 KSPConvergedReason kspreason; 142 143 PetscFunctionBegin; 144 snes->numFailures = 0; 145 snes->numLinearSolveFailures = 0; 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 /* set parameter for default relative tolerance convergence test */ 168 snes->ttol = fnorm*snes->rtol; 169 if (snes->ops->converged) { 170 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 171 } 172 if (snes->reason) PetscFunctionReturn(0); 173 174 for (i=0; i<maxits; i++) { 175 176 /* Call general purpose update function */ 177 if (snes->ops->update) { 178 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 179 } 180 181 /* Solve J Y = F, where J is Jacobian matrix */ 182 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 183 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 184 ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr); 185 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 186 if (kspreason < 0) { 187 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 188 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 189 PetscFunctionReturn(0); 190 } 191 } 192 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 193 snes->linear_its += lits; 194 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 195 196 if (neP->precheckstep) { 197 PetscTruth changed_y = PETSC_FALSE; 198 ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr); 199 } 200 201 if (PetscLogPrintInfo){ 202 ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 203 } 204 205 /* Compute a (scaled) negative update in the line search routine: 206 Y <- X - lambda*Y 207 and evaluate G = function(Y) (depends on the line search). 208 */ 209 ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 210 ynorm = 1; gnorm = fnorm; 211 ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 212 ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 213 TMP = F; F = G; snes->vec_func_always = F; G = TMP; 214 TMP = X; X = W; snes->vec_sol_always = X; W = TMP; 215 fnorm = gnorm; 216 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 217 218 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 219 snes->iter = i+1; 220 snes->norm = fnorm; 221 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 222 SNESLogConvHistory(snes,fnorm,lits); 223 SNESMonitor(snes,i+1,fnorm); 224 225 if (!lssucceed) { 226 PetscTruth ismin; 227 if (++snes->numFailures >= snes->maxFailures) { 228 snes->reason = SNES_DIVERGED_LS_FAILURE; 229 ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); 230 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 231 break; 232 } 233 } 234 235 /* Test for convergence */ 236 if (snes->ops->converged) { 237 ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 238 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 239 if (snes->reason) break; 240 } 241 } 242 if (X != snes->vec_sol) { 243 ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 244 } 245 if (F != snes->vec_func) { 246 ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 247 } 248 snes->vec_sol_always = snes->vec_sol; 249 snes->vec_func_always = snes->vec_func; 250 if (i == maxits) { 251 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 252 snes->reason = SNES_DIVERGED_MAX_IT; 253 } 254 PetscFunctionReturn(0); 255 } 256 /* -------------------------------------------------------------------------- */ 257 /* 258 SNESSetUp_LS - Sets up the internal data structures for the later use 259 of the SNESLS nonlinear solver. 260 261 Input Parameter: 262 . snes - the SNES context 263 . x - the solution vector 264 265 Application Interface Routine: SNESSetUp() 266 267 Notes: 268 For basic use of the SNES solvers, the user need not explicitly call 269 SNESSetUp(), since these actions will automatically occur during 270 the call to SNESSolve(). 271 */ 272 #undef __FUNCT__ 273 #define __FUNCT__ "SNESSetUp_LS" 274 PetscErrorCode SNESSetUp_LS(SNES snes) 275 { 276 PetscErrorCode ierr; 277 278 PetscFunctionBegin; 279 if (!snes->work) { 280 snes->nwork = 4; 281 ierr = VecDuplicateVecs(snes->vec_sol_always,snes->nwork,&snes->work);CHKERRQ(ierr); 282 ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 283 snes->vec_sol_update_always = snes->work[3]; 284 } 285 PetscFunctionReturn(0); 286 } 287 /* -------------------------------------------------------------------------- */ 288 /* 289 SNESDestroy_LS - Destroys the private SNES_LS context that was created 290 with SNESCreate_LS(). 291 292 Input Parameter: 293 . snes - the SNES context 294 295 Application Interface Routine: SNESDestroy() 296 */ 297 #undef __FUNCT__ 298 #define __FUNCT__ "SNESDestroy_LS" 299 PetscErrorCode SNESDestroy_LS(SNES snes) 300 { 301 PetscErrorCode ierr; 302 303 PetscFunctionBegin; 304 if (snes->nwork) { 305 ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 306 } 307 ierr = PetscFree(snes->data);CHKERRQ(ierr); 308 309 /* clear composed functions */ 310 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 311 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr); 312 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr); 313 314 PetscFunctionReturn(0); 315 } 316 /* -------------------------------------------------------------------------- */ 317 #undef __FUNCT__ 318 #define __FUNCT__ "SNESLineSearchNo" 319 320 /*@C 321 SNESLineSearchNo - This routine is not a line search at all; 322 it simply uses the full Newton step. Thus, this routine is intended 323 to serve as a template and is not recommended for general use. 324 325 Collective on SNES and Vec 326 327 Input Parameters: 328 + snes - nonlinear context 329 . lsctx - optional context for line search (not used here) 330 . x - current iterate 331 . f - residual evaluated at x 332 . y - search direction 333 . w - work vector 334 - fnorm - 2-norm of f 335 336 Output Parameters: 337 + g - residual evaluated at new iterate y 338 . w - new iterate 339 . gnorm - 2-norm of g 340 . ynorm - 2-norm of search length 341 - flag - PETSC_TRUE on success, PETSC_FALSE on failure 342 343 Options Database Key: 344 . -snes_ls basic - Activates SNESLineSearchNo() 345 346 Level: advanced 347 348 .keywords: SNES, nonlinear, line search, cubic 349 350 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 351 SNESLineSearchSet(), SNESLineSearchNoNorms() 352 @*/ 353 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) 354 { 355 PetscErrorCode ierr; 356 SNES_LS *neP = (SNES_LS*)snes->data; 357 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 358 359 PetscFunctionBegin; 360 *flag = PETSC_TRUE; 361 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 362 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 363 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 364 if (neP->postcheckstep) { 365 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 366 } 367 if (changed_y) { 368 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 369 } 370 ierr = SNESComputeFunction(snes,w,g); 371 if (PetscExceptionValue(ierr)) { 372 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 373 } 374 CHKERRQ(ierr); 375 376 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 377 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 378 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 379 PetscFunctionReturn(0); 380 } 381 /* -------------------------------------------------------------------------- */ 382 #undef __FUNCT__ 383 #define __FUNCT__ "SNESLineSearchNoNorms" 384 385 /*@C 386 SNESLineSearchNoNorms - This routine is not a line search at 387 all; it simply uses the full Newton step. This version does not 388 even compute the norm of the function or search direction; this 389 is intended only when you know the full step is fine and are 390 not checking for convergence of the nonlinear iteration (for 391 example, you are running always for a fixed number of Newton steps). 392 393 Collective on SNES and Vec 394 395 Input Parameters: 396 + snes - nonlinear context 397 . lsctx - optional context for line search (not used here) 398 . x - current iterate 399 . f - residual evaluated at x 400 . y - search direction 401 . w - work vector 402 - fnorm - 2-norm of f 403 404 Output Parameters: 405 + g - residual evaluated at new iterate y 406 . w - new iterate 407 . gnorm - not changed 408 . ynorm - not changed 409 - flag - set to PETSC_TRUE indicating a successful line search 410 411 Options Database Key: 412 . -snes_ls basicnonorms - Activates SNESLineSearchNoNorms() 413 414 Notes: 415 SNESLineSearchNoNorms() must be used in conjunction with 416 either the options 417 $ -snes_no_convergence_test -snes_max_it <its> 418 or alternatively a user-defined custom test set via 419 SNESSetConvergenceTest(); or a -snes_max_it of 1, 420 otherwise, the SNES solver will generate an error. 421 422 During the final iteration this will not evaluate the function at 423 the solution point. This is to save a function evaluation while 424 using pseudo-timestepping. 425 426 The residual norms printed by monitoring routines such as 427 SNESMonitorDefault() (as activated via -snes_monitor) will not be 428 correct, since they are not computed. 429 430 Level: advanced 431 432 .keywords: SNES, nonlinear, line search, cubic 433 434 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 435 SNESLineSearchSet(), SNESLineSearchNo() 436 @*/ 437 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) 438 { 439 PetscErrorCode ierr; 440 SNES_LS *neP = (SNES_LS*)snes->data; 441 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 442 443 PetscFunctionBegin; 444 *flag = PETSC_TRUE; 445 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 446 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 447 if (neP->postcheckstep) { 448 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 449 } 450 if (changed_y) { 451 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 452 } 453 454 /* don't evaluate function the last time through */ 455 if (snes->iter < snes->max_its-1) { 456 ierr = SNESComputeFunction(snes,w,g); 457 if (PetscExceptionValue(ierr)) { 458 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 459 } 460 CHKERRQ(ierr); 461 } 462 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 463 PetscFunctionReturn(0); 464 } 465 /* -------------------------------------------------------------------------- */ 466 #undef __FUNCT__ 467 #define __FUNCT__ "SNESLineSearchCubic" 468 /*@C 469 SNESLineSearchCubic - Performs a cubic line search (default line search method). 470 471 Collective on SNES 472 473 Input Parameters: 474 + snes - nonlinear context 475 . lsctx - optional context for line search (not used here) 476 . x - current iterate 477 . f - residual evaluated at x 478 . y - search direction 479 . w - work vector 480 - fnorm - 2-norm of f 481 482 Output Parameters: 483 + g - residual evaluated at new iterate y 484 . w - new iterate 485 . gnorm - 2-norm of g 486 . ynorm - 2-norm of search length 487 - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 488 489 Options Database Key: 490 $ -snes_ls cubic - Activates SNESLineSearchCubic() 491 492 Notes: 493 This line search is taken from "Numerical Methods for Unconstrained 494 Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 495 496 Level: advanced 497 498 .keywords: SNES, nonlinear, line search, cubic 499 500 .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 501 @*/ 502 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) 503 { 504 /* 505 Note that for line search purposes we work with with the related 506 minimization problem: 507 min z(x): R^n -> R, 508 where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 509 */ 510 511 PetscReal steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 512 PetscReal maxstep,minlambda,alpha,lambda,lambdatemp; 513 #if defined(PETSC_USE_COMPLEX) 514 PetscScalar cinitslope,clambda; 515 #endif 516 PetscErrorCode ierr; 517 PetscInt count; 518 SNES_LS *neP = (SNES_LS*)snes->data; 519 PetscScalar scale; 520 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 521 522 PetscFunctionBegin; 523 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 524 *flag = PETSC_TRUE; 525 alpha = neP->alpha; 526 maxstep = neP->maxstep; 527 steptol = neP->steptol; 528 529 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 530 if (!*ynorm) { 531 ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr); 532 *gnorm = fnorm; 533 ierr = VecCopy(x,w);CHKERRQ(ierr); 534 ierr = VecCopy(f,g);CHKERRQ(ierr); 535 *flag = PETSC_FALSE; 536 goto theend1; 537 } 538 if (*ynorm > maxstep) { /* Step too big, so scale back */ 539 scale = maxstep/(*ynorm); 540 ierr = PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",PetscRealPart(scale),*ynorm);CHKERRQ(ierr); 541 ierr = VecScale(y,scale);CHKERRQ(ierr); 542 *ynorm = maxstep; 543 } 544 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 545 minlambda = steptol/rellength; 546 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 547 #if defined(PETSC_USE_COMPLEX) 548 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 549 initslope = PetscRealPart(cinitslope); 550 #else 551 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 552 #endif 553 if (initslope > 0.0) initslope = -initslope; 554 if (initslope == 0.0) initslope = -1.0; 555 556 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 557 if (snes->nfuncs >= snes->max_funcs) { 558 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 559 *flag = PETSC_FALSE; 560 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 561 goto theend1; 562 } 563 ierr = SNESComputeFunction(snes,w,g); 564 if (PetscExceptionValue(ierr)) { 565 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 566 } 567 CHKERRQ(ierr); 568 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 569 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 570 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 571 ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr); 572 goto theend1; 573 } 574 575 /* Fit points with quadratic */ 576 lambda = 1.0; 577 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 578 lambdaprev = lambda; 579 gnormprev = *gnorm; 580 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 581 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 582 else lambda = lambdatemp; 583 584 #if defined(PETSC_USE_COMPLEX) 585 clambda = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr); 586 #else 587 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 588 #endif 589 if (snes->nfuncs >= snes->max_funcs) { 590 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 591 *flag = PETSC_FALSE; 592 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 593 goto theend1; 594 } 595 ierr = SNESComputeFunction(snes,w,g); 596 if (PetscExceptionValue(ierr)) { 597 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 598 } 599 CHKERRQ(ierr); 600 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 601 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 602 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 603 ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 604 goto theend1; 605 } 606 607 /* Fit points with cubic */ 608 count = 1; 609 while (count < 10000) { 610 if (lambda <= minlambda) { /* bad luck; use full step */ 611 ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr); 612 ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 613 *flag = PETSC_FALSE; 614 break; 615 } 616 t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 617 t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 618 a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 619 b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 620 d = b*b - 3*a*initslope; 621 if (d < 0.0) d = 0.0; 622 if (a == 0.0) { 623 lambdatemp = -initslope/(2.0*b); 624 } else { 625 lambdatemp = (-b + sqrt(d))/(3.0*a); 626 } 627 lambdaprev = lambda; 628 gnormprev = *gnorm; 629 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 630 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 631 else lambda = lambdatemp; 632 #if defined(PETSC_USE_COMPLEX) 633 clambda = lambda; 634 ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr); 635 #else 636 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 637 #endif 638 if (snes->nfuncs >= snes->max_funcs) { 639 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 640 ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 641 *flag = PETSC_FALSE; 642 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 643 break; 644 } 645 ierr = SNESComputeFunction(snes,w,g); 646 if (PetscExceptionValue(ierr)) { 647 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 648 } 649 CHKERRQ(ierr); 650 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 651 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 652 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */ 653 ierr = PetscInfo1(snes,"Cubically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 654 break; 655 } else { 656 ierr = PetscInfo1(snes,"Cubic step no good, shrinking lambda, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 657 } 658 count++; 659 } 660 if (count >= 10000) { 661 SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation"); 662 } 663 theend1: 664 /* Optional user-defined check for line search step validity */ 665 if (neP->postcheckstep && *flag) { 666 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 667 if (changed_y) { 668 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 669 } 670 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 671 ierr = SNESComputeFunction(snes,w,g); 672 if (PetscExceptionValue(ierr)) { 673 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 674 } 675 CHKERRQ(ierr); 676 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 677 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 678 ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr); 679 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 680 ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr); 681 } 682 } 683 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 684 PetscFunctionReturn(0); 685 } 686 /* -------------------------------------------------------------------------- */ 687 #undef __FUNCT__ 688 #define __FUNCT__ "SNESLineSearchQuadratic" 689 /*@C 690 SNESLineSearchQuadratic - Performs a quadratic line search. 691 692 Collective on SNES and Vec 693 694 Input Parameters: 695 + snes - the SNES context 696 . lsctx - optional context for line search (not used here) 697 . x - current iterate 698 . f - residual evaluated at x 699 . y - search direction 700 . w - work vector 701 - fnorm - 2-norm of f 702 703 Output Parameters: 704 + g - residual evaluated at new iterate w 705 . w - new iterate (x + alpha*y) 706 . gnorm - 2-norm of g 707 . ynorm - 2-norm of search length 708 - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 709 710 Options Database Keys: 711 + -snes_ls quadratic - Activates SNESLineSearchQuadratic() 712 . -snes_ls_alpha <alpha> - Sets alpha 713 . -snes_ls_maxstep <max> - Sets maxstep 714 - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 715 will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 716 the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 717 Notes: 718 Use SNESLineSearchSet() to set this routine within the SNESLS method. 719 720 Level: advanced 721 722 .keywords: SNES, nonlinear, quadratic, line search 723 724 .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 725 @*/ 726 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) 727 { 728 /* 729 Note that for line search purposes we work with with the related 730 minimization problem: 731 min z(x): R^n -> R, 732 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 733 */ 734 PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,rellength; 735 #if defined(PETSC_USE_COMPLEX) 736 PetscScalar cinitslope,clambda; 737 #endif 738 PetscErrorCode ierr; 739 PetscInt count; 740 SNES_LS *neP = (SNES_LS*)snes->data; 741 PetscScalar scale; 742 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 743 744 PetscFunctionBegin; 745 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 746 *flag = PETSC_TRUE; 747 alpha = neP->alpha; 748 maxstep = neP->maxstep; 749 steptol = neP->steptol; 750 751 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 752 if (*ynorm == 0.0) { 753 ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr); 754 *gnorm = fnorm; 755 ierr = VecCopy(x,w);CHKERRQ(ierr); 756 ierr = VecCopy(f,g);CHKERRQ(ierr); 757 *flag = PETSC_FALSE; 758 goto theend2; 759 } 760 if (*ynorm > maxstep) { /* Step too big, so scale back */ 761 scale = maxstep/(*ynorm); 762 ierr = VecScale(y,scale);CHKERRQ(ierr); 763 *ynorm = maxstep; 764 } 765 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 766 minlambda = steptol/rellength; 767 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 768 #if defined(PETSC_USE_COMPLEX) 769 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 770 initslope = PetscRealPart(cinitslope); 771 #else 772 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 773 #endif 774 if (initslope > 0.0) initslope = -initslope; 775 if (initslope == 0.0) initslope = -1.0; 776 777 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 778 if (snes->nfuncs >= snes->max_funcs) { 779 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 780 *flag = PETSC_FALSE; 781 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 782 goto theend2; 783 } 784 ierr = SNESComputeFunction(snes,w,g); 785 if (PetscExceptionValue(ierr)) { 786 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 787 } 788 CHKERRQ(ierr); 789 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 790 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 791 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 792 ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr); 793 goto theend2; 794 } 795 796 /* Fit points with quadratic */ 797 lambda = 1.0; 798 count = 1; 799 while (PETSC_TRUE) { 800 if (lambda <= minlambda) { /* bad luck; use full step */ 801 ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr); 802 ierr = PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 803 ierr = VecCopy(x,w);CHKERRQ(ierr); 804 *flag = PETSC_FALSE; 805 break; 806 } 807 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 808 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 809 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 810 else lambda = lambdatemp; 811 812 #if defined(PETSC_USE_COMPLEX) 813 clambda = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr); 814 #else 815 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 816 #endif 817 if (snes->nfuncs >= snes->max_funcs) { 818 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 819 ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 820 *flag = PETSC_FALSE; 821 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 822 break; 823 } 824 ierr = SNESComputeFunction(snes,w,g); 825 if (PetscExceptionValue(ierr)) { 826 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 827 } 828 CHKERRQ(ierr); 829 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 830 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 831 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 832 ierr = PetscInfo1(snes,"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(w,-1.0,y,x);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 before the line search is called. 1016 1017 Input Parameters: 1018 + snes - nonlinear context obtained from SNESCreate() 1019 . func - pointer to function 1020 - checkctx - optional user-defined context for use by step checking routine 1021 1022 Collective on SNES 1023 1024 Calling sequence of func: 1025 .vb 1026 int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y) 1027 .ve 1028 where func returns an error code of 0 on success and a nonzero 1029 on failure. 1030 1031 Input parameters for func: 1032 + snes - nonlinear context 1033 . checkctx - optional user-defined context for use by step checking routine 1034 . x - previous iterate 1035 - y - new search direction and length 1036 1037 Output parameters for func: 1038 + y - search direction (possibly changed) 1039 - changed_y - indicates search direction was changed by this routine 1040 1041 Level: advanced 1042 1043 Notes: All line searches accept the new iterate computed by the line search checking routine. 1044 1045 .keywords: SNES, nonlinear, set, line search check, step check, routine 1046 1047 .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 1048 @*/ 1049 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx) 1050 { 1051 PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*); 1052 1053 PetscFunctionBegin; 1054 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 1055 if (f) { 1056 ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 1057 } 1058 PetscFunctionReturn(0); 1059 } 1060 1061 /* -------------------------------------------------------------------------- */ 1062 typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/ 1063 typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*); /* force argument to next function to not be extern C*/ 1064 EXTERN_C_BEGIN 1065 #undef __FUNCT__ 1066 #define __FUNCT__ "SNESLineSearchSetPostCheck_LS" 1067 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx) 1068 { 1069 PetscFunctionBegin; 1070 ((SNES_LS *)(snes->data))->postcheckstep = func; 1071 ((SNES_LS *)(snes->data))->postcheck = checkctx; 1072 PetscFunctionReturn(0); 1073 } 1074 EXTERN_C_END 1075 1076 EXTERN_C_BEGIN 1077 #undef __FUNCT__ 1078 #define __FUNCT__ "SNESLineSearchSetPreCheck_LS" 1079 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx) 1080 { 1081 PetscFunctionBegin; 1082 ((SNES_LS *)(snes->data))->precheckstep = func; 1083 ((SNES_LS *)(snes->data))->precheck = checkctx; 1084 PetscFunctionReturn(0); 1085 } 1086 EXTERN_C_END 1087 1088 /* 1089 SNESView_LS - Prints info from the SNESLS data structure. 1090 1091 Input Parameters: 1092 . SNES - the SNES context 1093 . viewer - visualization context 1094 1095 Application Interface Routine: SNESView() 1096 */ 1097 #undef __FUNCT__ 1098 #define __FUNCT__ "SNESView_LS" 1099 static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 1100 { 1101 SNES_LS *ls = (SNES_LS *)snes->data; 1102 const char *cstr; 1103 PetscErrorCode ierr; 1104 PetscTruth iascii; 1105 1106 PetscFunctionBegin; 1107 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1108 if (iascii) { 1109 if (ls->LineSearch == SNESLineSearchNo) cstr = "SNESLineSearchNo"; 1110 else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic"; 1111 else if (ls->LineSearch == SNESLineSearchCubic) cstr = "SNESLineSearchCubic"; 1112 else cstr = "unknown"; 1113 ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1114 ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, steptol=%G\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 1115 } else { 1116 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 1117 } 1118 PetscFunctionReturn(0); 1119 } 1120 /* -------------------------------------------------------------------------- */ 1121 /* 1122 SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 1123 1124 Input Parameter: 1125 . snes - the SNES context 1126 1127 Application Interface Routine: SNESSetFromOptions() 1128 */ 1129 #undef __FUNCT__ 1130 #define __FUNCT__ "SNESSetFromOptions_LS" 1131 static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 1132 { 1133 SNES_LS *ls = (SNES_LS *)snes->data; 1134 const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1135 PetscErrorCode ierr; 1136 PetscInt indx; 1137 PetscTruth flg; 1138 1139 PetscFunctionBegin; 1140 ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 1141 ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 1142 ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 1143 ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr); 1144 1145 ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 1146 if (flg) { 1147 switch (indx) { 1148 case 0: 1149 ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 1150 break; 1151 case 1: 1152 ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1153 break; 1154 case 2: 1155 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr); 1156 break; 1157 case 3: 1158 ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr); 1159 break; 1160 } 1161 } 1162 ierr = PetscOptionsTail();CHKERRQ(ierr); 1163 PetscFunctionReturn(0); 1164 } 1165 /* -------------------------------------------------------------------------- */ 1166 #undef __FUNCT__ 1167 #define __FUNCT__ "SNESConverged_LS" 1168 /*@C 1169 SNESConverged_LS - Monitors the convergence of the line search 1170 method SNESLS for solving systems of nonlinear equations (default). 1171 1172 Collective on SNES 1173 1174 Input Parameters: 1175 + snes - the SNES context 1176 . it - the iteration (0 indicates before any Newton steps) 1177 . xnorm - 2-norm of current iterate 1178 . pnorm - 2-norm of current step 1179 . fnorm - 2-norm of function at current iterate 1180 - dummy - unused context 1181 1182 Output Parameter: 1183 . reason - one of 1184 $ SNES_CONVERGED_FNORM_ABS - (fnorm < abstol), 1185 $ SNES_CONVERGED_PNORM_RELATIVE - (pnorm < xtol*xnorm), 1186 $ SNES_CONVERGED_FNORM_RELATIVE - (fnorm < rtol*fnorm0), 1187 $ SNES_DIVERGED_FUNCTION_COUNT - (nfct > maxf), 1188 $ SNES_DIVERGED_FNORM_NAN - (fnorm == NaN), 1189 $ SNES_CONVERGED_ITERATING - (otherwise), 1190 1191 where 1192 + maxf - maximum number of function evaluations, 1193 set with SNESSetTolerances() 1194 . nfct - number of function evaluations, 1195 . abstol - absolute function norm tolerance, 1196 set with SNESSetTolerances() 1197 - rtol - relative function norm tolerance, set with SNESSetTolerances() 1198 1199 Level: intermediate 1200 1201 .keywords: SNES, nonlinear, default, converged, convergence 1202 1203 .seealso: SNESSetConvergenceTest() 1204 @*/ 1205 PetscErrorCode PETSCSNES_DLLEXPORT SNESConverged_LS(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 1206 { 1207 PetscErrorCode ierr; 1208 PetscFunctionBegin; 1209 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1210 PetscValidType(snes,1); 1211 PetscValidPointer(reason,6); 1212 ierr = SNESDefaultConverged(snes,it,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr); 1213 PetscFunctionReturn(0); 1214 } 1215 /* -------------------------------------------------------------------------- */ 1216 /*MC 1217 SNESLS - Newton based nonlinear solver that uses a line search 1218 1219 Options Database: 1220 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1221 . -snes_ls_alpha <alpha> - Sets alpha 1222 . -snes_ls_maxstep <max> - Sets maxstep 1223 - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 1224 will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 1225 the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 1226 1227 Notes: This is the default nonlinear solver in SNES 1228 1229 Level: beginner 1230 1231 .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 1232 SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 1233 SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck() 1234 1235 M*/ 1236 EXTERN_C_BEGIN 1237 #undef __FUNCT__ 1238 #define __FUNCT__ "SNESCreate_LS" 1239 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes) 1240 { 1241 PetscErrorCode ierr; 1242 SNES_LS *neP; 1243 1244 PetscFunctionBegin; 1245 snes->ops->setup = SNESSetUp_LS; 1246 snes->ops->solve = SNESSolve_LS; 1247 snes->ops->destroy = SNESDestroy_LS; 1248 snes->ops->converged = SNESConverged_LS; 1249 snes->ops->setfromoptions = SNESSetFromOptions_LS; 1250 snes->ops->view = SNESView_LS; 1251 snes->nwork = 0; 1252 1253 ierr = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr); 1254 snes->data = (void*)neP; 1255 neP->alpha = 1.e-4; 1256 neP->maxstep = 1.e8; 1257 neP->steptol = 1.e-12; 1258 neP->LineSearch = SNESLineSearchCubic; 1259 neP->lsP = PETSC_NULL; 1260 neP->postcheckstep = PETSC_NULL; 1261 neP->postcheck = PETSC_NULL; 1262 neP->precheckstep = PETSC_NULL; 1263 neP->precheck = PETSC_NULL; 1264 1265 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C", 1266 "SNESLineSearchSet_LS", 1267 SNESLineSearchSet_LS);CHKERRQ(ierr); 1268 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C", 1269 "SNESLineSearchSetPostCheck_LS", 1270 SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr); 1271 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C", 1272 "SNESLineSearchSetPreCheck_LS", 1273 SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr); 1274 1275 PetscFunctionReturn(0); 1276 } 1277 EXTERN_C_END 1278 1279 1280 1281