1 2 #include <../src/snes/impls/ls/lsimpl.h> 3 4 /* 5 Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function, 6 || F(u) ||_2 but not a zero, F(u) = 0. In the case when one cannot compute J^T F we use the fact that 7 0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More 8 for this trick. One assumes that the probability that W is in the null space of J is very, very small. 9 */ 10 #undef __FUNCT__ 11 #define __FUNCT__ "SNESLSCheckLocalMin_Private" 12 PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool *ismin) 13 { 14 PetscReal a1; 15 PetscErrorCode ierr; 16 PetscBool hastranspose; 17 18 PetscFunctionBegin; 19 *ismin = PETSC_FALSE; 20 ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 21 if (hastranspose) { 22 /* Compute || J^T F|| */ 23 ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr); 24 ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr); 25 ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr); 26 if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE; 27 } else { 28 Vec work; 29 PetscScalar result; 30 PetscReal wnorm; 31 32 ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr); 33 ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr); 34 ierr = VecDuplicate(W,&work);CHKERRQ(ierr); 35 ierr = MatMult(A,W,work);CHKERRQ(ierr); 36 ierr = VecDot(F,work,&result);CHKERRQ(ierr); 37 ierr = VecDestroy(work);CHKERRQ(ierr); 38 a1 = PetscAbsScalar(result)/(fnorm*wnorm); 39 ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr); 40 if (a1 < 1.e-4) *ismin = PETSC_TRUE; 41 } 42 PetscFunctionReturn(0); 43 } 44 45 /* 46 Checks if J^T(F - J*X) = 0 47 */ 48 #undef __FUNCT__ 49 #define __FUNCT__ "SNESLSCheckResidual_Private" 50 PetscErrorCode SNESLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2) 51 { 52 PetscReal a1,a2; 53 PetscErrorCode ierr; 54 PetscBool hastranspose; 55 56 PetscFunctionBegin; 57 ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr); 58 if (hastranspose) { 59 ierr = MatMult(A,X,W1);CHKERRQ(ierr); 60 ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr); 61 62 /* Compute || J^T W|| */ 63 ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr); 64 ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr); 65 ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr); 66 if (a1 != 0.0) { 67 ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr); 68 } 69 } 70 PetscFunctionReturn(0); 71 } 72 73 /* -------------------------------------------------------------------- 74 75 This file implements a truncated Newton method with a line search, 76 for solving a system of nonlinear equations, using the KSP, Vec, 77 and Mat interfaces for linear solvers, vectors, and matrices, 78 respectively. 79 80 The following basic routines are required for each nonlinear solver: 81 SNESCreate_XXX() - Creates a nonlinear solver context 82 SNESSetFromOptions_XXX() - Sets runtime options 83 SNESSolve_XXX() - Solves the nonlinear system 84 SNESDestroy_XXX() - Destroys the nonlinear solver context 85 The suffix "_XXX" denotes a particular implementation, in this case 86 we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving 87 systems of nonlinear equations with a line search (LS) method. 88 These routines are actually called via the common user interface 89 routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and 90 SNESDestroy(), so the application code interface remains identical 91 for all nonlinear solvers. 92 93 Another key routine is: 94 SNESSetUp_XXX() - Prepares for the use of a nonlinear solver 95 by setting data structures and options. The interface routine SNESSetUp() 96 is not usually called directly by the user, but instead is called by 97 SNESSolve() if necessary. 98 99 Additional basic routines are: 100 SNESView_XXX() - Prints details of runtime options that 101 have actually been used. 102 These are called by application codes via the interface routines 103 SNESView(). 104 105 The various types of solvers (preconditioners, Krylov subspace methods, 106 nonlinear solvers, timesteppers) are all organized similarly, so the 107 above description applies to these categories also. 108 109 -------------------------------------------------------------------- */ 110 /* 111 SNESSolve_LS - Solves a nonlinear system with a truncated Newton 112 method with a line search. 113 114 Input Parameters: 115 . snes - the SNES context 116 117 Output Parameter: 118 . outits - number of iterations until termination 119 120 Application Interface Routine: SNESSolve() 121 122 Notes: 123 This implements essentially a truncated Newton method with a 124 line search. By default a cubic backtracking line search 125 is employed, as described in the text "Numerical Methods for 126 Unconstrained Optimization and Nonlinear Equations" by Dennis 127 and Schnabel. 128 */ 129 #undef __FUNCT__ 130 #define __FUNCT__ "SNESSolve_LS" 131 PetscErrorCode SNESSolve_LS(SNES snes) 132 { 133 SNES_LS *neP = (SNES_LS*)snes->data; 134 PetscErrorCode ierr; 135 PetscInt maxits,i,lits; 136 PetscBool lssucceed; 137 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 138 PetscReal fnorm,gnorm,xnorm=0,ynorm; 139 Vec Y,X,F,G,W; 140 KSPConvergedReason kspreason; 141 142 PetscFunctionBegin; 143 snes->numFailures = 0; 144 snes->numLinearSolveFailures = 0; 145 snes->reason = SNES_CONVERGED_ITERATING; 146 147 maxits = snes->max_its; /* maximum number of iterations */ 148 X = snes->vec_sol; /* solution vector */ 149 F = snes->vec_func; /* residual vector */ 150 Y = snes->work[0]; /* work vectors */ 151 G = snes->work[1]; 152 W = snes->work[2]; 153 154 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 155 snes->iter = 0; 156 snes->norm = 0.0; 157 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 158 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 159 if (snes->domainerror) { 160 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 161 PetscFunctionReturn(0); 162 } 163 ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 164 ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- ||x|| */ 165 ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr); 166 ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr); 167 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 168 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 169 snes->norm = fnorm; 170 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 171 SNESLogConvHistory(snes,fnorm,0); 172 SNESMonitor(snes,0,fnorm); 173 174 /* set parameter for default relative tolerance convergence test */ 175 snes->ttol = fnorm*snes->rtol; 176 /* test convergence */ 177 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 178 if (snes->reason) PetscFunctionReturn(0); 179 180 for (i=0; i<maxits; i++) { 181 182 /* Call general purpose update function */ 183 if (snes->ops->update) { 184 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 185 } 186 187 /* Solve J Y = F, where J is Jacobian matrix */ 188 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 189 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr); 190 ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr); 191 ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); 192 if (kspreason < 0) { 193 if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { 194 ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); 195 snes->reason = SNES_DIVERGED_LINEAR_SOLVE; 196 break; 197 } 198 } 199 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 200 snes->linear_its += lits; 201 ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr); 202 203 if (neP->precheckstep) { 204 PetscBool changed_y = PETSC_FALSE; 205 ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr); 206 } 207 208 if (PetscLogPrintInfo){ 209 ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr); 210 } 211 212 /* Compute a (scaled) negative update in the line search routine: 213 Y <- X - lambda*Y 214 and evaluate G = function(Y) (depends on the line search). 215 */ 216 ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 217 ynorm = 1; gnorm = fnorm; 218 ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 219 ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr); 220 if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break; 221 if (snes->domainerror) { 222 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 223 PetscFunctionReturn(0); 224 } 225 if (!lssucceed) { 226 if (++snes->numFailures >= snes->maxFailures) { 227 PetscBool ismin; 228 snes->reason = SNES_DIVERGED_LINE_SEARCH; 229 ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr); 230 if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN; 231 break; 232 } 233 } 234 /* Update function and solution vectors */ 235 fnorm = gnorm; 236 ierr = VecCopy(G,F);CHKERRQ(ierr); 237 ierr = VecCopy(W,X);CHKERRQ(ierr); 238 /* Monitor convergence */ 239 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 240 snes->iter = i+1; 241 snes->norm = fnorm; 242 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 243 SNESLogConvHistory(snes,snes->norm,lits); 244 SNESMonitor(snes,snes->iter,snes->norm); 245 /* Test for convergence, xnorm = || X || */ 246 if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); } 247 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 248 if (snes->reason) break; 249 } 250 if (i == maxits) { 251 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr); 252 if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 253 } 254 PetscFunctionReturn(0); 255 } 256 /* -------------------------------------------------------------------------- */ 257 /* 258 SNESSetUp_LS - Sets up the internal data structures for the later use 259 of the SNESLS nonlinear solver. 260 261 Input Parameter: 262 . snes - the SNES context 263 . x - the solution vector 264 265 Application Interface Routine: SNESSetUp() 266 267 Notes: 268 For basic use of the SNES solvers, the user need not explicitly call 269 SNESSetUp(), since these actions will automatically occur during 270 the call to SNESSolve(). 271 */ 272 #undef __FUNCT__ 273 #define __FUNCT__ "SNESSetUp_LS" 274 PetscErrorCode SNESSetUp_LS(SNES snes) 275 { 276 PetscErrorCode ierr; 277 278 PetscFunctionBegin; 279 ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr); 280 PetscFunctionReturn(0); 281 } 282 /* -------------------------------------------------------------------------- */ 283 /* 284 SNESDestroy_LS - Destroys the private SNES_LS context that was created 285 with SNESCreate_LS(). 286 287 Input Parameter: 288 . snes - the SNES context 289 290 Application Interface Routine: SNESDestroy() 291 */ 292 #undef __FUNCT__ 293 #define __FUNCT__ "SNESDestroy_LS" 294 PetscErrorCode SNESDestroy_LS(SNES snes) 295 { 296 SNES_LS *ls = (SNES_LS*) snes->data; 297 PetscErrorCode ierr; 298 299 PetscFunctionBegin; 300 if (ls->monitor) { 301 ierr = PetscViewerASCIIMonitorDestroy(ls->monitor);CHKERRQ(ierr); 302 } 303 ierr = PetscFree(snes->data);CHKERRQ(ierr); 304 305 /* clear composed functions */ 306 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr); 307 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetParams_C","",PETSC_NULL);CHKERRQ(ierr); 308 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 309 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr); 310 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr); 311 312 PetscFunctionReturn(0); 313 } 314 /* -------------------------------------------------------------------------- */ 315 #undef __FUNCT__ 316 #define __FUNCT__ "SNESLineSearchNo" 317 318 /*@C 319 SNESLineSearchNo - This routine is not a line search at all; 320 it simply uses the full Newton step. Thus, this routine is intended 321 to serve as a template and is not recommended for general use. 322 323 Logically Collective on SNES and Vec 324 325 Input Parameters: 326 + snes - nonlinear context 327 . lsctx - optional context for line search (not used here) 328 . x - current iterate 329 . f - residual evaluated at x 330 . y - search direction 331 . fnorm - 2-norm of f 332 - xnorm - norm of x if known, otherwise 0 333 334 Output Parameters: 335 + g - residual evaluated at new iterate y 336 . w - new iterate 337 . gnorm - 2-norm of g 338 . ynorm - 2-norm of search length 339 - flag - PETSC_TRUE on success, PETSC_FALSE on failure 340 341 Options Database Key: 342 . -snes_ls basic - Activates SNESLineSearchNo() 343 344 Level: advanced 345 346 .keywords: SNES, nonlinear, line search, cubic 347 348 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 349 SNESLineSearchSet(), SNESLineSearchNoNorms() 350 @*/ 351 PetscErrorCode SNESLineSearchNo(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 352 { 353 PetscErrorCode ierr; 354 SNES_LS *neP = (SNES_LS*)snes->data; 355 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 356 357 PetscFunctionBegin; 358 *flag = PETSC_TRUE; 359 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 360 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 361 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 362 if (neP->postcheckstep) { 363 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 364 } 365 if (changed_y) { 366 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 367 } 368 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 369 if (!snes->domainerror) { 370 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 371 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 372 } 373 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 374 PetscFunctionReturn(0); 375 } 376 /* -------------------------------------------------------------------------- */ 377 #undef __FUNCT__ 378 #define __FUNCT__ "SNESLineSearchNoNorms" 379 380 /*@C 381 SNESLineSearchNoNorms - This routine is not a line search at 382 all; it simply uses the full Newton step. This version does not 383 even compute the norm of the function or search direction; this 384 is intended only when you know the full step is fine and are 385 not checking for convergence of the nonlinear iteration (for 386 example, you are running always for a fixed number of Newton steps). 387 388 Logically Collective on SNES and Vec 389 390 Input Parameters: 391 + snes - nonlinear context 392 . lsctx - optional context for line search (not used here) 393 . x - current iterate 394 . f - residual evaluated at x 395 . y - search direction 396 . w - work vector 397 . fnorm - 2-norm of f 398 - xnorm - norm of x if known, otherwise 0 399 400 Output Parameters: 401 + g - residual evaluated at new iterate y 402 . w - new iterate 403 . gnorm - not changed 404 . ynorm - not changed 405 - flag - set to PETSC_TRUE indicating a successful line search 406 407 Options Database Key: 408 . -snes_ls basicnonorms - Activates SNESLineSearchNoNorms() 409 410 Notes: 411 SNESLineSearchNoNorms() must be used in conjunction with 412 either the options 413 $ -snes_no_convergence_test -snes_max_it <its> 414 or alternatively a user-defined custom test set via 415 SNESSetConvergenceTest(); or a -snes_max_it of 1, 416 otherwise, the SNES solver will generate an error. 417 418 During the final iteration this will not evaluate the function at 419 the solution point. This is to save a function evaluation while 420 using pseudo-timestepping. 421 422 The residual norms printed by monitoring routines such as 423 SNESMonitorDefault() (as activated via -snes_monitor) will not be 424 correct, since they are not computed. 425 426 Level: advanced 427 428 .keywords: SNES, nonlinear, line search, cubic 429 430 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 431 SNESLineSearchSet(), SNESLineSearchNo() 432 @*/ 433 PetscErrorCode SNESLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 434 { 435 PetscErrorCode ierr; 436 SNES_LS *neP = (SNES_LS*)snes->data; 437 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 438 439 PetscFunctionBegin; 440 *flag = PETSC_TRUE; 441 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 442 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 443 if (neP->postcheckstep) { 444 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 445 } 446 if (changed_y) { 447 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 448 } 449 450 /* don't evaluate function the last time through */ 451 if (snes->iter < snes->max_its-1) { 452 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 453 } 454 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 455 PetscFunctionReturn(0); 456 } 457 /* -------------------------------------------------------------------------- */ 458 #undef __FUNCT__ 459 #define __FUNCT__ "SNESLineSearchCubic" 460 /*@C 461 SNESLineSearchCubic - Performs a cubic line search (default line search method). 462 463 Collective on SNES 464 465 Input Parameters: 466 + snes - nonlinear context 467 . lsctx - optional context for line search (not used here) 468 . x - current iterate 469 . f - residual evaluated at x 470 . y - search direction 471 . w - work vector 472 . fnorm - 2-norm of f 473 - xnorm - norm of x if known, otherwise 0 474 475 Output Parameters: 476 + g - residual evaluated at new iterate y 477 . w - new iterate 478 . gnorm - 2-norm of g 479 . ynorm - 2-norm of search length 480 - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 481 482 Options Database Key: 483 + -snes_ls cubic - Activates SNESLineSearchCubic() 484 . -snes_ls_alpha <alpha> - Sets alpha 485 . -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y) 486 - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 487 488 489 Notes: 490 This line search is taken from "Numerical Methods for Unconstrained 491 Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 492 493 Level: advanced 494 495 .keywords: SNES, nonlinear, line search, cubic 496 497 .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 498 @*/ 499 PetscErrorCode SNESLineSearchCubic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 500 { 501 /* 502 Note that for line search purposes we work with with the related 503 minimization problem: 504 min z(x): R^n -> R, 505 where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 506 */ 507 508 PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 509 PetscReal minlambda,lambda,lambdatemp; 510 #if defined(PETSC_USE_COMPLEX) 511 PetscScalar cinitslope; 512 #endif 513 PetscErrorCode ierr; 514 PetscInt count; 515 SNES_LS *neP = (SNES_LS*)snes->data; 516 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 517 MPI_Comm comm; 518 519 PetscFunctionBegin; 520 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 521 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 522 *flag = PETSC_TRUE; 523 524 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 525 if (*ynorm == 0.0) { 526 if (neP->monitor) { 527 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 528 } 529 *gnorm = fnorm; 530 ierr = VecCopy(x,w);CHKERRQ(ierr); 531 ierr = VecCopy(f,g);CHKERRQ(ierr); 532 *flag = PETSC_FALSE; 533 goto theend1; 534 } 535 if (*ynorm > neP->maxstep) { /* Step too big, so scale back */ 536 if (neP->monitor) { 537 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Scaling step by %G old ynorm %G\n",neP->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 538 } 539 ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr); 540 *ynorm = neP->maxstep; 541 } 542 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 543 minlambda = neP->minlambda/rellength; 544 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 545 #if defined(PETSC_USE_COMPLEX) 546 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 547 initslope = PetscRealPart(cinitslope); 548 #else 549 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 550 #endif 551 if (initslope > 0.0) initslope = -initslope; 552 if (initslope == 0.0) initslope = -1.0; 553 554 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 555 if (snes->nfuncs >= snes->max_funcs) { 556 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 557 *flag = PETSC_FALSE; 558 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 559 goto theend1; 560 } 561 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 562 if (snes->domainerror) { 563 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 564 PetscFunctionReturn(0); 565 } 566 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 567 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 568 ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 569 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */ 570 if (neP->monitor) { 571 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 572 } 573 goto theend1; 574 } 575 576 /* Fit points with quadratic */ 577 lambda = 1.0; 578 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 579 lambdaprev = lambda; 580 gnormprev = *gnorm; 581 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 582 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 583 else lambda = lambdatemp; 584 585 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 586 if (snes->nfuncs >= snes->max_funcs) { 587 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 588 *flag = PETSC_FALSE; 589 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 590 goto theend1; 591 } 592 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 593 if (snes->domainerror) { 594 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 595 PetscFunctionReturn(0); 596 } 597 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 598 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 599 if (neP->monitor) { 600 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 601 } 602 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */ 603 if (neP->monitor) { 604 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 605 } 606 goto theend1; 607 } 608 609 /* Fit points with cubic */ 610 count = 1; 611 while (PETSC_TRUE) { 612 if (lambda <= minlambda) { 613 if (neP->monitor) { 614 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 615 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)minlambda,(double)lambda,(double)initslope);CHKERRQ(ierr); 616 } 617 *flag = PETSC_FALSE; 618 break; 619 } 620 t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 621 t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 622 a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 623 b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 624 d = b*b - 3*a*initslope; 625 if (d < 0.0) d = 0.0; 626 if (a == 0.0) { 627 lambdatemp = -initslope/(2.0*b); 628 } else { 629 lambdatemp = (-b + sqrt(d))/(3.0*a); 630 } 631 lambdaprev = lambda; 632 gnormprev = *gnorm; 633 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 634 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 635 else lambda = lambdatemp; 636 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 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);CHKERRQ(ierr); 645 if (snes->domainerror) { 646 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 647 PetscFunctionReturn(0); 648 } 649 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 650 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 651 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* is reduction enough? */ 652 if (neP->monitor) { 653 ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,(double)lambda);CHKERRQ(ierr); 654 } 655 break; 656 } else { 657 if (neP->monitor) { 658 ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,(double)lambda);CHKERRQ(ierr); 659 } 660 } 661 count++; 662 } 663 theend1: 664 /* Optional user-defined check for line search step validity */ 665 if (neP->postcheckstep && *flag) { 666 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 667 if (changed_y) { 668 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 669 } 670 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 671 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 672 if (snes->domainerror) { 673 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 674 PetscFunctionReturn(0); 675 } 676 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 677 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 678 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 679 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 680 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 681 } 682 } 683 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 684 PetscFunctionReturn(0); 685 } 686 /* -------------------------------------------------------------------------- */ 687 #undef __FUNCT__ 688 #define __FUNCT__ "SNESLineSearchQuadratic" 689 /*@C 690 SNESLineSearchQuadratic - Performs a quadratic line search. 691 692 Collective on SNES and Vec 693 694 Input Parameters: 695 + snes - the SNES context 696 . lsctx - optional context for line search (not used here) 697 . x - current iterate 698 . f - residual evaluated at x 699 . y - search direction 700 . w - work vector 701 . fnorm - 2-norm of f 702 - xnorm - norm of x if known, otherwise 0 703 704 Output Parameters: 705 + g - residual evaluated at new iterate w 706 . w - new iterate (x + lambda*y) 707 . gnorm - 2-norm of g 708 . ynorm - 2-norm of search length 709 - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 710 711 Options Database Keys: 712 + -snes_ls quadratic - Activates SNESLineSearchQuadratic() 713 . -snes_ls_alpha <alpha> - Sets alpha 714 . -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y) 715 - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 716 717 Notes: 718 Use SNESLineSearchSet() to set this routine within the SNESLS method. 719 720 Level: advanced 721 722 .keywords: SNES, nonlinear, quadratic, line search 723 724 .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 725 @*/ 726 PetscErrorCode SNESLineSearchQuadratic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 727 { 728 /* 729 Note that for line search purposes we work with with the related 730 minimization problem: 731 min z(x): R^n -> R, 732 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 733 */ 734 PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 735 #if defined(PETSC_USE_COMPLEX) 736 PetscScalar cinitslope; 737 #endif 738 PetscErrorCode ierr; 739 PetscInt count; 740 SNES_LS *neP = (SNES_LS*)snes->data; 741 PetscBool 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 747 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 748 if (*ynorm == 0.0) { 749 if (neP->monitor) { 750 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 751 } 752 *gnorm = fnorm; 753 ierr = VecCopy(x,w);CHKERRQ(ierr); 754 ierr = VecCopy(f,g);CHKERRQ(ierr); 755 *flag = PETSC_FALSE; 756 goto theend2; 757 } 758 if (*ynorm > neP->maxstep) { /* Step too big, so scale back */ 759 ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr); 760 *ynorm = neP->maxstep; 761 } 762 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 763 minlambda = neP->minlambda/rellength; 764 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 765 #if defined(PETSC_USE_COMPLEX) 766 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 767 initslope = PetscRealPart(cinitslope); 768 #else 769 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 770 #endif 771 if (initslope > 0.0) initslope = -initslope; 772 if (initslope == 0.0) initslope = -1.0; 773 ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 774 775 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 776 if (snes->nfuncs >= snes->max_funcs) { 777 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 778 *flag = PETSC_FALSE; 779 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 780 goto theend2; 781 } 782 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 783 if (snes->domainerror) { 784 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 785 PetscFunctionReturn(0); 786 } 787 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 788 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 789 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */ 790 if (neP->monitor) { 791 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line Search: Using full step\n");CHKERRQ(ierr); 792 } 793 goto theend2; 794 } 795 796 /* Fit points with quadratic */ 797 lambda = 1.0; 798 count = 1; 799 while (PETSC_TRUE) { 800 if (lambda <= minlambda) { /* bad luck; use full step */ 801 if (neP->monitor) { 802 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 803 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"Line search: fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr); 804 } 805 ierr = VecCopy(x,w);CHKERRQ(ierr); 806 *flag = PETSC_FALSE; 807 break; 808 } 809 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 810 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 811 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 812 else lambda = lambdatemp; 813 814 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 815 if (snes->nfuncs >= snes->max_funcs) { 816 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 817 ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr); 818 *flag = PETSC_FALSE; 819 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 820 break; 821 } 822 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 823 if (snes->domainerror) { 824 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 825 PetscFunctionReturn(0); 826 } 827 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 828 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 829 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */ 830 if (neP->monitor) { 831 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 832 } 833 break; 834 } 835 count++; 836 } 837 theend2: 838 /* Optional user-defined check for line search step validity */ 839 if (neP->postcheckstep) { 840 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 841 if (changed_y) { 842 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 843 } 844 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 845 ierr = SNESComputeFunction(snes,w,g); 846 if (snes->domainerror) { 847 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 848 PetscFunctionReturn(0); 849 } 850 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 851 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 852 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 853 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 854 if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 855 } 856 } 857 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 858 PetscFunctionReturn(0); 859 } 860 861 /* -------------------------------------------------------------------------- */ 862 #undef __FUNCT__ 863 #define __FUNCT__ "SNESLineSearchSet" 864 /*@C 865 SNESLineSearchSet - Sets the line search routine to be used 866 by the method SNESLS. 867 868 Input Parameters: 869 + snes - nonlinear context obtained from SNESCreate() 870 . lsctx - optional user-defined context for use by line search 871 - func - pointer to int function 872 873 Logically Collective on SNES 874 875 Available Routines: 876 + SNESLineSearchCubic() - default line search 877 . SNESLineSearchQuadratic() - quadratic line search 878 . SNESLineSearchNo() - the full Newton step (actually not a line search) 879 - SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 880 881 Options Database Keys: 882 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 883 . -snes_ls_alpha <alpha> - Sets alpha 884 . -snes_ls_maxstep <maxstep> - Sets maximum step the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y) 885 - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 886 887 Calling sequence of func: 888 .vb 889 func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool *flag) 890 .ve 891 892 Input parameters for func: 893 + snes - nonlinear context 894 . lsctx - optional user-defined context for line search 895 . x - current iterate 896 . f - residual evaluated at x 897 . y - search direction 898 - fnorm - 2-norm of f 899 900 Output parameters for func: 901 + g - residual evaluated at new iterate y 902 . w - new iterate 903 . gnorm - 2-norm of g 904 . ynorm - 2-norm of search length 905 - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 906 907 Level: advanced 908 909 .keywords: SNES, nonlinear, set, line search, routine 910 911 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(), 912 SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck() 913 @*/ 914 PetscErrorCode SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *),void *lsctx) 915 { 916 PetscErrorCode ierr; 917 918 PetscFunctionBegin; 919 ierr = PetscTryMethod(snes,"SNESLineSearchSet_C",(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *),void*),(snes,func,lsctx));CHKERRQ(ierr); 920 PetscFunctionReturn(0); 921 } 922 923 typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/ 924 /* -------------------------------------------------------------------------- */ 925 EXTERN_C_BEGIN 926 #undef __FUNCT__ 927 #define __FUNCT__ "SNESLineSearchSet_LS" 928 PetscErrorCode SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx) 929 { 930 PetscFunctionBegin; 931 ((SNES_LS *)(snes->data))->LineSearch = func; 932 ((SNES_LS *)(snes->data))->lsP = lsctx; 933 PetscFunctionReturn(0); 934 } 935 EXTERN_C_END 936 /* -------------------------------------------------------------------------- */ 937 #undef __FUNCT__ 938 #define __FUNCT__ "SNESLineSearchSetPostCheck" 939 /*@C 940 SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed 941 by the line search routine in the Newton-based method SNESLS. 942 943 Input Parameters: 944 + snes - nonlinear context obtained from SNESCreate() 945 . func - pointer to function 946 - checkctx - optional user-defined context for use by step checking routine 947 948 Logically Collective on SNES 949 950 Calling sequence of func: 951 .vb 952 int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscBool *changed_y,PetscBool *changed_w) 953 .ve 954 where func returns an error code of 0 on success and a nonzero 955 on failure. 956 957 Input parameters for func: 958 + snes - nonlinear context 959 . checkctx - optional user-defined context for use by step checking routine 960 . x - previous iterate 961 . y - new search direction and length 962 - w - current candidate iterate 963 964 Output parameters for func: 965 + y - search direction (possibly changed) 966 . w - current iterate (possibly modified) 967 . changed_y - indicates search direction was changed by this routine 968 - changed_w - indicates current iterate was changed by this routine 969 970 Level: advanced 971 972 Notes: All line searches accept the new iterate computed by the line search checking routine. 973 974 Only one of changed_y and changed_w can be PETSC_TRUE 975 976 On input w = x + y 977 978 SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control 979 to the checking routine, and then (3) compute the corresponding nonlinear 980 function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 981 982 SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a 983 candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 984 routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 985 were made to the candidate iterate in the checking routine (as indicated 986 by flag=PETSC_TRUE). The overhead of this extra function re-evaluation can be 987 very costly, so use this feature with caution! 988 989 .keywords: SNES, nonlinear, set, line search check, step check, routine 990 991 .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck() 992 @*/ 993 PetscErrorCode SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void *checkctx) 994 { 995 PetscErrorCode ierr; 996 997 PetscFunctionBegin; 998 ierr = PetscTryMethod(snes,"SNESLineSearchSetPostCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void*),(snes,func,checkctx));CHKERRQ(ierr); 999 PetscFunctionReturn(0); 1000 } 1001 1002 #undef __FUNCT__ 1003 #define __FUNCT__ "SNESLineSearchSetPreCheck" 1004 /*@C 1005 SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve 1006 before the line search is called. 1007 1008 Input Parameters: 1009 + snes - nonlinear context obtained from SNESCreate() 1010 . func - pointer to function 1011 - checkctx - optional user-defined context for use by step checking routine 1012 1013 Logically Collective on SNES 1014 1015 Calling sequence of func: 1016 .vb 1017 int func (SNES snes, Vec x,Vec y,void *checkctx, PetscBool *changed_y) 1018 .ve 1019 where func returns an error code of 0 on success and a nonzero 1020 on failure. 1021 1022 Input parameters for func: 1023 + snes - nonlinear context 1024 . checkctx - optional user-defined context for use by step checking routine 1025 . x - previous iterate 1026 - y - new search direction and length 1027 1028 Output parameters for func: 1029 + y - search direction (possibly changed) 1030 - changed_y - indicates search direction was changed by this routine 1031 1032 Level: advanced 1033 1034 Notes: All line searches accept the new iterate computed by the line search checking routine. 1035 1036 .keywords: SNES, nonlinear, set, line search check, step check, routine 1037 1038 .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 1039 @*/ 1040 PetscErrorCode SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscBool *),void *checkctx) 1041 { 1042 PetscErrorCode ierr; 1043 1044 PetscFunctionBegin; 1045 ierr = PetscTryMethod(snes,"SNESLineSearchSetPreCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscBool *),void*),(snes,func,checkctx));CHKERRQ(ierr); 1046 PetscFunctionReturn(0); 1047 } 1048 1049 #undef __FUNCT__ 1050 #define __FUNCT__ "SNESLineSearchSetMonitor" 1051 /*@C 1052 SNESLineSearchSetMonitor - Prints information about the progress or lack of progress of the line search 1053 1054 Input Parameters: 1055 + snes - nonlinear context obtained from SNESCreate() 1056 - flg - PETSC_TRUE to monitor the line search 1057 1058 Logically Collective on SNES 1059 1060 Options Database: 1061 . -snes_ls_monitor 1062 1063 Level: intermediate 1064 1065 1066 .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 1067 @*/ 1068 PetscErrorCode SNESLineSearchSetMonitor(SNES snes,PetscBool flg) 1069 { 1070 PetscErrorCode ierr; 1071 1072 PetscFunctionBegin; 1073 ierr = PetscTryMethod(snes,"SNESLineSearchSetMonitor_C",(SNES,PetscBool),(snes,flg));CHKERRQ(ierr); 1074 PetscFunctionReturn(0); 1075 } 1076 1077 /* -------------------------------------------------------------------------- */ 1078 typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *); /* force argument to next function to not be extern C*/ 1079 typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscBool *); /* force argument to next function to not be extern C*/ 1080 EXTERN_C_BEGIN 1081 #undef __FUNCT__ 1082 #define __FUNCT__ "SNESLineSearchSetPostCheck_LS" 1083 PetscErrorCode SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx) 1084 { 1085 PetscFunctionBegin; 1086 ((SNES_LS *)(snes->data))->postcheckstep = func; 1087 ((SNES_LS *)(snes->data))->postcheck = checkctx; 1088 PetscFunctionReturn(0); 1089 } 1090 EXTERN_C_END 1091 1092 EXTERN_C_BEGIN 1093 #undef __FUNCT__ 1094 #define __FUNCT__ "SNESLineSearchSetPreCheck_LS" 1095 PetscErrorCode SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx) 1096 { 1097 PetscFunctionBegin; 1098 ((SNES_LS *)(snes->data))->precheckstep = func; 1099 ((SNES_LS *)(snes->data))->precheck = checkctx; 1100 PetscFunctionReturn(0); 1101 } 1102 EXTERN_C_END 1103 1104 EXTERN_C_BEGIN 1105 #undef __FUNCT__ 1106 #define __FUNCT__ "SNESLineSearchSetMonitor_LS" 1107 PetscErrorCode SNESLineSearchSetMonitor_LS(SNES snes,PetscBool flg) 1108 { 1109 SNES_LS *ls = (SNES_LS*)snes->data; 1110 PetscErrorCode ierr; 1111 1112 PetscFunctionBegin; 1113 if (flg && !ls->monitor) { 1114 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&ls->monitor);CHKERRQ(ierr); 1115 } else if (!flg && ls->monitor) { 1116 ierr = PetscViewerASCIIMonitorDestroy(ls->monitor);CHKERRQ(ierr); 1117 } 1118 PetscFunctionReturn(0); 1119 } 1120 EXTERN_C_END 1121 1122 /* 1123 SNESView_LS - Prints info from the SNESLS data structure. 1124 1125 Input Parameters: 1126 . SNES - the SNES context 1127 . viewer - visualization context 1128 1129 Application Interface Routine: SNESView() 1130 */ 1131 #undef __FUNCT__ 1132 #define __FUNCT__ "SNESView_LS" 1133 static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 1134 { 1135 SNES_LS *ls = (SNES_LS *)snes->data; 1136 const char *cstr; 1137 PetscErrorCode ierr; 1138 PetscBool iascii; 1139 1140 PetscFunctionBegin; 1141 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1142 if (iascii) { 1143 if (ls->LineSearch == SNESLineSearchNo) cstr = "SNESLineSearchNo"; 1144 else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic"; 1145 else if (ls->LineSearch == SNESLineSearchCubic) cstr = "SNESLineSearchCubic"; 1146 else cstr = "unknown"; 1147 ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1148 ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",ls->alpha,ls->maxstep,ls->minlambda);CHKERRQ(ierr); 1149 } else { 1150 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 1151 } 1152 PetscFunctionReturn(0); 1153 } 1154 /* -------------------------------------------------------------------------- */ 1155 /* 1156 SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 1157 1158 Input Parameter: 1159 . snes - the SNES context 1160 1161 Application Interface Routine: SNESSetFromOptions() 1162 */ 1163 #undef __FUNCT__ 1164 #define __FUNCT__ "SNESSetFromOptions_LS" 1165 static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 1166 { 1167 SNES_LS *ls = (SNES_LS *)snes->data; 1168 const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1169 PetscErrorCode ierr; 1170 PetscInt indx; 1171 PetscBool flg,set; 1172 1173 PetscFunctionBegin; 1174 ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 1175 ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 1176 ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 1177 ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",ls->minlambda,&ls->minlambda,0);CHKERRQ(ierr); 1178 ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",ls->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 1179 if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 1180 1181 ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 1182 if (flg) { 1183 switch (indx) { 1184 case 0: 1185 ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 1186 break; 1187 case 1: 1188 ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1189 break; 1190 case 2: 1191 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr); 1192 break; 1193 case 3: 1194 ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr); 1195 break; 1196 } 1197 } 1198 ierr = PetscOptionsTail();CHKERRQ(ierr); 1199 PetscFunctionReturn(0); 1200 } 1201 1202 EXTERN_C_BEGIN 1203 extern PetscErrorCode SNESLineSearchSetParams_LS(SNES,PetscReal,PetscReal,PetscReal); 1204 EXTERN_C_END 1205 1206 /* -------------------------------------------------------------------------- */ 1207 /*MC 1208 SNESLS - Newton based nonlinear solver that uses a line search 1209 1210 Options Database: 1211 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1212 . -snes_ls_alpha <alpha> - Sets alpha 1213 . -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y) 1214 . -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 1215 - -snes_ls_monitor - print information about progress of line searches 1216 1217 1218 Notes: This is the default nonlinear solver in SNES 1219 1220 Level: beginner 1221 1222 .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 1223 SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 1224 SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 1225 1226 M*/ 1227 EXTERN_C_BEGIN 1228 #undef __FUNCT__ 1229 #define __FUNCT__ "SNESCreate_LS" 1230 PetscErrorCode SNESCreate_LS(SNES snes) 1231 { 1232 PetscErrorCode ierr; 1233 SNES_LS *neP; 1234 1235 PetscFunctionBegin; 1236 snes->ops->setup = SNESSetUp_LS; 1237 snes->ops->solve = SNESSolve_LS; 1238 snes->ops->destroy = SNESDestroy_LS; 1239 snes->ops->setfromoptions = SNESSetFromOptions_LS; 1240 snes->ops->view = SNESView_LS; 1241 snes->ops->reset = 0; 1242 1243 ierr = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr); 1244 snes->data = (void*)neP; 1245 neP->alpha = 1.e-4; 1246 neP->maxstep = 1.e8; 1247 neP->minlambda = 1.e-12; 1248 neP->LineSearch = SNESLineSearchCubic; 1249 neP->lsP = PETSC_NULL; 1250 neP->postcheckstep = PETSC_NULL; 1251 neP->postcheck = PETSC_NULL; 1252 neP->precheckstep = PETSC_NULL; 1253 neP->precheck = PETSC_NULL; 1254 1255 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_LS",SNESLineSearchSetMonitor_LS);CHKERRQ(ierr); 1256 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetParams_C","SNESLineSearchSetParams_LS",SNESLineSearchSetParams_LS);CHKERRQ(ierr); 1257 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);CHKERRQ(ierr); 1258 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr); 1259 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr); 1260 1261 PetscFunctionReturn(0); 1262 } 1263 EXTERN_C_END 1264 1265 1266 1267