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 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 167 if (snes->reason) PetscFunctionReturn(0); 168 169 for (i=0; i<maxits; i++) { 170 171 /* Call general purpose update function */ 172 if (snes->ops->update) { 173 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 174 } 175 176 /* Solve J Y = F, where J is Jacobian matrix */ 177 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 178 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 179 ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr); 180 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 181 if (kspreason < 0) { 182 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 183 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 184 PetscFunctionReturn(0); 185 } 186 } 187 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 188 189 if (neP->precheckstep) { 190 PetscTruth changed_y = PETSC_FALSE; 191 ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr); 192 } 193 194 if (PetscLogPrintInfo){ 195 ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 196 } 197 198 /* should check what happened to the linear solve? */ 199 snes->linear_its += lits; 200 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 201 202 /* Compute a (scaled) negative update in the line search routine: 203 Y <- X - lambda*Y 204 and evaluate G = function(Y) (depends on the line search). 205 */ 206 ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr); 207 ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 208 ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 209 TMP = F; F = G; snes->vec_func_always = F; G = TMP; 210 TMP = X; X = W; snes->vec_sol_always = X; W = TMP; 211 fnorm = gnorm; 212 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 213 214 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 215 snes->iter = i+1; 216 snes->norm = fnorm; 217 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 218 SNESLogConvHistory(snes,fnorm,lits); 219 SNESMonitor(snes,i+1,fnorm); 220 221 if (!lssucceed) { 222 PetscTruth ismin; 223 224 if (++snes->numFailures >= snes->maxFailures) { 225 snes->reason = SNES_DIVERGED_LS_FAILURE; 226 ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr); 227 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 228 break; 229 } 230 } 231 232 /* Test for convergence */ 233 if (snes->ops->converged) { 234 ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm = || X || */ 235 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 236 if (snes->reason) { 237 break; 238 } 239 } 240 } 241 if (X != snes->vec_sol) { 242 ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr); 243 } 244 if (F != snes->vec_func) { 245 ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr); 246 } 247 snes->vec_sol_always = snes->vec_sol; 248 snes->vec_func_always = snes->vec_func; 249 if (i == maxits) { 250 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 251 snes->reason = SNES_DIVERGED_MAX_IT; 252 } 253 PetscFunctionReturn(0); 254 } 255 /* -------------------------------------------------------------------------- */ 256 /* 257 SNESSetUp_LS - Sets up the internal data structures for the later use 258 of the SNESLS nonlinear solver. 259 260 Input Parameter: 261 . snes - the SNES context 262 . x - the solution vector 263 264 Application Interface Routine: SNESSetUp() 265 266 Notes: 267 For basic use of the SNES solvers, the user need not explicitly call 268 SNESSetUp(), since these actions will automatically occur during 269 the call to SNESSolve(). 270 */ 271 #undef __FUNCT__ 272 #define __FUNCT__ "SNESSetUp_LS" 273 PetscErrorCode SNESSetUp_LS(SNES snes) 274 { 275 PetscErrorCode ierr; 276 277 PetscFunctionBegin; 278 if (!snes->work) { 279 snes->nwork = 4; 280 ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr); 281 ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr); 282 snes->vec_sol_update_always = snes->work[3]; 283 } 284 PetscFunctionReturn(0); 285 } 286 /* -------------------------------------------------------------------------- */ 287 /* 288 SNESDestroy_LS - Destroys the private SNES_LS context that was created 289 with SNESCreate_LS(). 290 291 Input Parameter: 292 . snes - the SNES context 293 294 Application Interface Routine: SNESDestroy() 295 */ 296 #undef __FUNCT__ 297 #define __FUNCT__ "SNESDestroy_LS" 298 PetscErrorCode SNESDestroy_LS(SNES snes) 299 { 300 PetscErrorCode ierr; 301 302 PetscFunctionBegin; 303 if (snes->nwork) { 304 ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr); 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 Key: 710 . -snes_ls quadratic - Activates SNESLineSearchQuadratic() 711 712 Notes: 713 Use SNESLineSearchSet() to set this routine within the SNESLS method. 714 715 Level: advanced 716 717 .keywords: SNES, nonlinear, quadratic, line search 718 719 .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 720 @*/ 721 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) 722 { 723 /* 724 Note that for line search purposes we work with with the related 725 minimization problem: 726 min z(x): R^n -> R, 727 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 728 */ 729 PetscReal steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,rellength; 730 #if defined(PETSC_USE_COMPLEX) 731 PetscScalar cinitslope,clambda; 732 #endif 733 PetscErrorCode ierr; 734 PetscInt count; 735 SNES_LS *neP = (SNES_LS*)snes->data; 736 PetscScalar scale; 737 PetscTruth changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 738 739 PetscFunctionBegin; 740 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 741 *flag = PETSC_TRUE; 742 alpha = neP->alpha; 743 maxstep = neP->maxstep; 744 steptol = neP->steptol; 745 746 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 747 if (*ynorm == 0.0) { 748 ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr); 749 *gnorm = fnorm; 750 ierr = VecCopy(x,w);CHKERRQ(ierr); 751 ierr = VecCopy(f,g);CHKERRQ(ierr); 752 *flag = PETSC_FALSE; 753 goto theend2; 754 } 755 if (*ynorm > maxstep) { /* Step too big, so scale back */ 756 scale = maxstep/(*ynorm); 757 ierr = VecScale(y,scale);CHKERRQ(ierr); 758 *ynorm = maxstep; 759 } 760 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 761 minlambda = steptol/rellength; 762 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 763 #if defined(PETSC_USE_COMPLEX) 764 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 765 initslope = PetscRealPart(cinitslope); 766 #else 767 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 768 #endif 769 if (initslope > 0.0) initslope = -initslope; 770 if (initslope == 0.0) initslope = -1.0; 771 772 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 773 if (snes->nfuncs >= snes->max_funcs) { 774 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 775 *flag = PETSC_FALSE; 776 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 777 goto theend2; 778 } 779 ierr = SNESComputeFunction(snes,w,g); 780 if (PetscExceptionValue(ierr)) { 781 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 782 } 783 CHKERRQ(ierr); 784 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 785 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 786 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */ 787 ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr); 788 goto theend2; 789 } 790 791 /* Fit points with quadratic */ 792 lambda = 1.0; 793 count = 1; 794 while (PETSC_TRUE) { 795 if (lambda <= minlambda) { /* bad luck; use full step */ 796 ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr); 797 ierr = PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 798 ierr = VecCopy(x,w);CHKERRQ(ierr); 799 *flag = PETSC_FALSE; 800 break; 801 } 802 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 803 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 804 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 805 else lambda = lambdatemp; 806 807 #if defined(PETSC_USE_COMPLEX) 808 clambda = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr); 809 #else 810 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 811 #endif 812 if (snes->nfuncs >= snes->max_funcs) { 813 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 814 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); 815 *flag = PETSC_FALSE; 816 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 817 break; 818 } 819 ierr = SNESComputeFunction(snes,w,g); 820 if (PetscExceptionValue(ierr)) { 821 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 822 } 823 CHKERRQ(ierr); 824 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 825 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 826 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */ 827 ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 828 break; 829 } 830 count++; 831 } 832 theend2: 833 /* Optional user-defined check for line search step validity */ 834 if (neP->postcheckstep) { 835 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 836 if (changed_y) { 837 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 838 } 839 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 840 ierr = SNESComputeFunction(snes,w,g); 841 if (PetscExceptionValue(ierr)) { 842 PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr); 843 } 844 CHKERRQ(ierr); 845 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 846 ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr); 847 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 848 ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr); 849 if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 850 } 851 } 852 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 853 PetscFunctionReturn(0); 854 } 855 856 /* -------------------------------------------------------------------------- */ 857 #undef __FUNCT__ 858 #define __FUNCT__ "SNESLineSearchSet" 859 /*@C 860 SNESLineSearchSet - Sets the line search routine to be used 861 by the method SNESLS. 862 863 Input Parameters: 864 + snes - nonlinear context obtained from SNESCreate() 865 . lsctx - optional user-defined context for use by line search 866 - func - pointer to int function 867 868 Collective on SNES 869 870 Available Routines: 871 + SNESLineSearchCubic() - default line search 872 . SNESLineSearchQuadratic() - quadratic line search 873 . SNESLineSearchNo() - the full Newton step (actually not a line search) 874 - SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 875 876 Options Database Keys: 877 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 878 . -snes_ls_alpha <alpha> - Sets alpha 879 . -snes_ls_maxstep <max> - Sets maxstep 880 - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 881 will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 882 the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 883 884 Calling sequence of func: 885 .vb 886 func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w, 887 PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag) 888 .ve 889 890 Input parameters for func: 891 + snes - nonlinear context 892 . lsctx - optional user-defined context for line search 893 . x - current iterate 894 . f - residual evaluated at x 895 . y - search direction 896 . w - work vector 897 - fnorm - 2-norm of f 898 899 Output parameters for func: 900 + g - residual evaluated at new iterate y 901 . w - new iterate 902 . gnorm - 2-norm of g 903 . ynorm - 2-norm of search length 904 - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 905 906 Level: advanced 907 908 .keywords: SNES, nonlinear, set, line search, routine 909 910 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(), 911 SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck() 912 @*/ 913 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx) 914 { 915 PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*); 916 917 PetscFunctionBegin; 918 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr); 919 if (f) { 920 ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr); 921 } 922 PetscFunctionReturn(0); 923 } 924 925 typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/ 926 /* -------------------------------------------------------------------------- */ 927 EXTERN_C_BEGIN 928 #undef __FUNCT__ 929 #define __FUNCT__ "SNESLineSearchSet_LS" 930 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx) 931 { 932 PetscFunctionBegin; 933 ((SNES_LS *)(snes->data))->LineSearch = func; 934 ((SNES_LS *)(snes->data))->lsP = lsctx; 935 PetscFunctionReturn(0); 936 } 937 EXTERN_C_END 938 /* -------------------------------------------------------------------------- */ 939 #undef __FUNCT__ 940 #define __FUNCT__ "SNESLineSearchSetPostCheck" 941 /*@C 942 SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed 943 by the line search routine in the Newton-based method SNESLS. 944 945 Input Parameters: 946 + snes - nonlinear context obtained from SNESCreate() 947 . func - pointer to function 948 - checkctx - optional user-defined context for use by step checking routine 949 950 Collective on SNES 951 952 Calling sequence of func: 953 .vb 954 int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w) 955 .ve 956 where func returns an error code of 0 on success and a nonzero 957 on failure. 958 959 Input parameters for func: 960 + snes - nonlinear context 961 . checkctx - optional user-defined context for use by step checking routine 962 . x - previous iterate 963 . y - new search direction and length 964 - w - current candidate iterate 965 966 Output parameters for func: 967 + y - search direction (possibly changed) 968 . w - current iterate (possibly modified) 969 . changed_y - indicates search direction was changed by this routine 970 - changed_w - indicates current iterate was changed by this routine 971 972 Level: advanced 973 974 Notes: All line searches accept the new iterate computed by the line search checking routine. 975 976 Only one of changed_y and changed_w can be PETSC_TRUE 977 978 On input w = x + y 979 980 SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control 981 to the checking routine, and then (3) compute the corresponding nonlinear 982 function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 983 984 SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a 985 candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 986 routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 987 were made to the candidate iterate in the checking routine (as indicated 988 by flag=PETSC_TRUE). The overhead of this extra function re-evaluation can be 989 very costly, so use this feature with caution! 990 991 .keywords: SNES, nonlinear, set, line search check, step check, routine 992 993 .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck() 994 @*/ 995 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx) 996 { 997 PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*); 998 999 PetscFunctionBegin; 1000 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 1001 if (f) { 1002 ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 1003 } 1004 PetscFunctionReturn(0); 1005 } 1006 #undef __FUNCT__ 1007 #define __FUNCT__ "SNESLineSearchSetPreCheck" 1008 /*@C 1009 SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve 1010 before the line search is called. 1011 1012 Input Parameters: 1013 + snes - nonlinear context obtained from SNESCreate() 1014 . func - pointer to function 1015 - checkctx - optional user-defined context for use by step checking routine 1016 1017 Collective on SNES 1018 1019 Calling sequence of func: 1020 .vb 1021 int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y) 1022 .ve 1023 where func returns an error code of 0 on success and a nonzero 1024 on failure. 1025 1026 Input parameters for func: 1027 + snes - nonlinear context 1028 . checkctx - optional user-defined context for use by step checking routine 1029 . x - previous iterate 1030 - y - new search direction and length 1031 1032 Output parameters for func: 1033 + y - search direction (possibly changed) 1034 - changed_y - indicates search direction was changed by this routine 1035 1036 Level: advanced 1037 1038 Notes: All line searches accept the new iterate computed by the line search checking routine. 1039 1040 .keywords: SNES, nonlinear, set, line search check, step check, routine 1041 1042 .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 1043 @*/ 1044 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx) 1045 { 1046 PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*); 1047 1048 PetscFunctionBegin; 1049 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr); 1050 if (f) { 1051 ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr); 1052 } 1053 PetscFunctionReturn(0); 1054 } 1055 1056 /* -------------------------------------------------------------------------- */ 1057 typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/ 1058 typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*); /* force argument to next function to not be extern C*/ 1059 EXTERN_C_BEGIN 1060 #undef __FUNCT__ 1061 #define __FUNCT__ "SNESLineSearchSetPostCheck_LS" 1062 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx) 1063 { 1064 PetscFunctionBegin; 1065 ((SNES_LS *)(snes->data))->postcheckstep = func; 1066 ((SNES_LS *)(snes->data))->postcheck = checkctx; 1067 PetscFunctionReturn(0); 1068 } 1069 EXTERN_C_END 1070 1071 EXTERN_C_BEGIN 1072 #undef __FUNCT__ 1073 #define __FUNCT__ "SNESLineSearchSetPreCheck_LS" 1074 PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx) 1075 { 1076 PetscFunctionBegin; 1077 ((SNES_LS *)(snes->data))->precheckstep = func; 1078 ((SNES_LS *)(snes->data))->precheck = checkctx; 1079 PetscFunctionReturn(0); 1080 } 1081 EXTERN_C_END 1082 1083 /* 1084 SNESView_LS - Prints info from the SNESLS data structure. 1085 1086 Input Parameters: 1087 . SNES - the SNES context 1088 . viewer - visualization context 1089 1090 Application Interface Routine: SNESView() 1091 */ 1092 #undef __FUNCT__ 1093 #define __FUNCT__ "SNESView_LS" 1094 static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 1095 { 1096 SNES_LS *ls = (SNES_LS *)snes->data; 1097 const char *cstr; 1098 PetscErrorCode ierr; 1099 PetscTruth iascii; 1100 1101 PetscFunctionBegin; 1102 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1103 if (iascii) { 1104 if (ls->LineSearch == SNESLineSearchNo) cstr = "SNESLineSearchNo"; 1105 else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic"; 1106 else if (ls->LineSearch == SNESLineSearchCubic) cstr = "SNESLineSearchCubic"; 1107 else cstr = "unknown"; 1108 ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1109 ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, steptol=%G\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr); 1110 } else { 1111 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 1112 } 1113 PetscFunctionReturn(0); 1114 } 1115 /* -------------------------------------------------------------------------- */ 1116 /* 1117 SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 1118 1119 Input Parameter: 1120 . snes - the SNES context 1121 1122 Application Interface Routine: SNESSetFromOptions() 1123 */ 1124 #undef __FUNCT__ 1125 #define __FUNCT__ "SNESSetFromOptions_LS" 1126 static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 1127 { 1128 SNES_LS *ls = (SNES_LS *)snes->data; 1129 const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1130 PetscErrorCode ierr; 1131 PetscInt indx; 1132 PetscTruth flg; 1133 1134 PetscFunctionBegin; 1135 ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 1136 ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 1137 ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 1138 ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr); 1139 1140 ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 1141 if (flg) { 1142 switch (indx) { 1143 case 0: 1144 ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 1145 break; 1146 case 1: 1147 ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1148 break; 1149 case 2: 1150 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr); 1151 break; 1152 case 3: 1153 ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr); 1154 break; 1155 } 1156 } 1157 ierr = PetscOptionsTail();CHKERRQ(ierr); 1158 PetscFunctionReturn(0); 1159 } 1160 /* -------------------------------------------------------------------------- */ 1161 /*MC 1162 SNESLS - Newton based nonlinear solver that uses a line search 1163 1164 Options Database: 1165 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1166 . -snes_ls_alpha <alpha> - Sets alpha 1167 . -snes_ls_maxstep <max> - Sets maxstep 1168 - -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code 1169 will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length 1170 the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x) 1171 1172 Notes: This is the default nonlinear solver in SNES 1173 1174 Level: beginner 1175 1176 .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 1177 SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 1178 SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck() 1179 1180 M*/ 1181 EXTERN_C_BEGIN 1182 #undef __FUNCT__ 1183 #define __FUNCT__ "SNESCreate_LS" 1184 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes) 1185 { 1186 PetscErrorCode ierr; 1187 SNES_LS *neP; 1188 1189 PetscFunctionBegin; 1190 snes->ops->setup = SNESSetUp_LS; 1191 snes->ops->solve = SNESSolve_LS; 1192 snes->ops->destroy = SNESDestroy_LS; 1193 snes->ops->converged = SNESConverged_LS; 1194 snes->ops->setfromoptions = SNESSetFromOptions_LS; 1195 snes->ops->view = SNESView_LS; 1196 snes->nwork = 0; 1197 1198 ierr = PetscNew(SNES_LS,&neP);CHKERRQ(ierr); 1199 ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr); 1200 snes->data = (void*)neP; 1201 neP->alpha = 1.e-4; 1202 neP->maxstep = 1.e8; 1203 neP->steptol = 1.e-12; 1204 neP->LineSearch = SNESLineSearchCubic; 1205 neP->lsP = PETSC_NULL; 1206 neP->postcheckstep = PETSC_NULL; 1207 neP->postcheck = PETSC_NULL; 1208 neP->precheckstep = PETSC_NULL; 1209 neP->precheck = PETSC_NULL; 1210 1211 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C", 1212 "SNESLineSearchSet_LS", 1213 SNESLineSearchSet_LS);CHKERRQ(ierr); 1214 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C", 1215 "SNESLineSearchSetPostCheck_LS", 1216 SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr); 1217 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C", 1218 "SNESLineSearchSetPreCheck_LS", 1219 SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr); 1220 1221 PetscFunctionReturn(0); 1222 } 1223 EXTERN_C_END 1224 1225 1226 1227