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