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