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