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