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