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->work) { 314 ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr); 315 } 316 if (ls->monitor) { 317 ierr = PetscViewerASCIIMonitorDestroy(ls->monitor);CHKERRQ(ierr); 318 } 319 ierr = PetscFree(snes->data);CHKERRQ(ierr); 320 321 /* clear composed functions */ 322 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr); 323 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr); 324 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr); 325 326 PetscFunctionReturn(0); 327 } 328 /* -------------------------------------------------------------------------- */ 329 #undef __FUNCT__ 330 #define __FUNCT__ "SNESLineSearchNo" 331 332 /*@C 333 SNESLineSearchNo - This routine is not a line search at all; 334 it simply uses the full Newton step. Thus, this routine is intended 335 to serve as a template and is not recommended for general use. 336 337 Logically Collective on SNES and Vec 338 339 Input Parameters: 340 + snes - nonlinear context 341 . lsctx - optional context for line search (not used here) 342 . x - current iterate 343 . f - residual evaluated at x 344 . y - search direction 345 . fnorm - 2-norm of f 346 - xnorm - norm of x if known, otherwise 0 347 348 Output Parameters: 349 + g - residual evaluated at new iterate y 350 . w - new iterate 351 . gnorm - 2-norm of g 352 . ynorm - 2-norm of search length 353 - flag - PETSC_TRUE on success, PETSC_FALSE on failure 354 355 Options Database Key: 356 . -snes_ls basic - Activates SNESLineSearchNo() 357 358 Level: advanced 359 360 .keywords: SNES, nonlinear, line search, cubic 361 362 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 363 SNESLineSearchSet(), SNESLineSearchNoNorms() 364 @*/ 365 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) 366 { 367 PetscErrorCode ierr; 368 SNES_LS *neP = (SNES_LS*)snes->data; 369 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 370 371 PetscFunctionBegin; 372 *flag = PETSC_TRUE; 373 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 374 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); /* ynorm = || y || */ 375 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 376 if (neP->postcheckstep) { 377 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 378 } 379 if (changed_y) { 380 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 381 } 382 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 383 if (!snes->domainerror) { 384 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); /* gnorm = || g || */ 385 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 386 } 387 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 388 PetscFunctionReturn(0); 389 } 390 /* -------------------------------------------------------------------------- */ 391 #undef __FUNCT__ 392 #define __FUNCT__ "SNESLineSearchNoNorms" 393 394 /*@C 395 SNESLineSearchNoNorms - This routine is not a line search at 396 all; it simply uses the full Newton step. This version does not 397 even compute the norm of the function or search direction; this 398 is intended only when you know the full step is fine and are 399 not checking for convergence of the nonlinear iteration (for 400 example, you are running always for a fixed number of Newton steps). 401 402 Logically Collective on SNES and Vec 403 404 Input Parameters: 405 + snes - nonlinear context 406 . lsctx - optional context for line search (not used here) 407 . x - current iterate 408 . f - residual evaluated at x 409 . y - search direction 410 . w - work vector 411 . fnorm - 2-norm of f 412 - xnorm - norm of x if known, otherwise 0 413 414 Output Parameters: 415 + g - residual evaluated at new iterate y 416 . w - new iterate 417 . gnorm - not changed 418 . ynorm - not changed 419 - flag - set to PETSC_TRUE indicating a successful line search 420 421 Options Database Key: 422 . -snes_ls basicnonorms - Activates SNESLineSearchNoNorms() 423 424 Notes: 425 SNESLineSearchNoNorms() must be used in conjunction with 426 either the options 427 $ -snes_no_convergence_test -snes_max_it <its> 428 or alternatively a user-defined custom test set via 429 SNESSetConvergenceTest(); or a -snes_max_it of 1, 430 otherwise, the SNES solver will generate an error. 431 432 During the final iteration this will not evaluate the function at 433 the solution point. This is to save a function evaluation while 434 using pseudo-timestepping. 435 436 The residual norms printed by monitoring routines such as 437 SNESMonitorDefault() (as activated via -snes_monitor) will not be 438 correct, since they are not computed. 439 440 Level: advanced 441 442 .keywords: SNES, nonlinear, line search, cubic 443 444 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), 445 SNESLineSearchSet(), SNESLineSearchNo() 446 @*/ 447 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) 448 { 449 PetscErrorCode ierr; 450 SNES_LS *neP = (SNES_LS*)snes->data; 451 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 452 453 PetscFunctionBegin; 454 *flag = PETSC_TRUE; 455 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 456 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 457 if (neP->postcheckstep) { 458 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 459 } 460 if (changed_y) { 461 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); /* w <- x - y */ 462 } 463 464 /* don't evaluate function the last time through */ 465 if (snes->iter < snes->max_its-1) { 466 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 467 } 468 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 469 PetscFunctionReturn(0); 470 } 471 /* -------------------------------------------------------------------------- */ 472 #undef __FUNCT__ 473 #define __FUNCT__ "SNESLineSearchCubic" 474 /*@C 475 SNESLineSearchCubic - Performs a cubic line search (default line search method). 476 477 Collective on SNES 478 479 Input Parameters: 480 + snes - nonlinear context 481 . lsctx - optional context for line search (not used here) 482 . x - current iterate 483 . f - residual evaluated at x 484 . y - search direction 485 . w - work vector 486 . fnorm - 2-norm of f 487 - xnorm - norm of x if known, otherwise 0 488 489 Output Parameters: 490 + g - residual evaluated at new iterate y 491 . w - new iterate 492 . gnorm - 2-norm of g 493 . ynorm - 2-norm of search length 494 - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 495 496 Options Database Key: 497 + -snes_ls cubic - Activates SNESLineSearchCubic() 498 . -snes_ls_alpha <alpha> - Sets alpha 499 . -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) 500 - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 501 502 503 Notes: 504 This line search is taken from "Numerical Methods for Unconstrained 505 Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325. 506 507 Level: advanced 508 509 .keywords: SNES, nonlinear, line search, cubic 510 511 .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 512 @*/ 513 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) 514 { 515 /* 516 Note that for line search purposes we work with with the related 517 minimization problem: 518 min z(x): R^n -> R, 519 where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2. 520 */ 521 522 PetscReal initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength; 523 PetscReal minlambda,lambda,lambdatemp; 524 #if defined(PETSC_USE_COMPLEX) 525 PetscScalar cinitslope; 526 #endif 527 PetscErrorCode ierr; 528 PetscInt count; 529 SNES_LS *neP = (SNES_LS*)snes->data; 530 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 531 MPI_Comm comm; 532 533 PetscFunctionBegin; 534 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 535 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 536 *flag = PETSC_TRUE; 537 538 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 539 if (*ynorm == 0.0) { 540 if (neP->monitor) { 541 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); 542 } 543 *gnorm = fnorm; 544 ierr = VecCopy(x,w);CHKERRQ(ierr); 545 ierr = VecCopy(f,g);CHKERRQ(ierr); 546 *flag = PETSC_FALSE; 547 goto theend1; 548 } 549 if (*ynorm > neP->maxstep) { /* Step too big, so scale back */ 550 if (neP->monitor) { 551 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Scaling step by %G old ynorm %G\n",neP->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr); 552 } 553 ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr); 554 *ynorm = neP->maxstep; 555 } 556 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 557 minlambda = neP->minlambda/rellength; 558 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 559 #if defined(PETSC_USE_COMPLEX) 560 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 561 initslope = PetscRealPart(cinitslope); 562 #else 563 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 564 #endif 565 if (initslope > 0.0) initslope = -initslope; 566 if (initslope == 0.0) initslope = -1.0; 567 568 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 569 if (snes->nfuncs >= snes->max_funcs) { 570 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 571 *flag = PETSC_FALSE; 572 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 573 goto theend1; 574 } 575 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 576 if (snes->domainerror) { 577 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 578 PetscFunctionReturn(0); 579 } 580 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 581 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 582 ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 583 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */ 584 if (neP->monitor) { 585 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr); 586 } 587 goto theend1; 588 } 589 590 /* Fit points with quadratic */ 591 lambda = 1.0; 592 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 593 lambdaprev = lambda; 594 gnormprev = *gnorm; 595 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 596 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 597 else lambda = lambdatemp; 598 599 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 600 if (snes->nfuncs >= snes->max_funcs) { 601 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); 602 *flag = PETSC_FALSE; 603 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 604 goto theend1; 605 } 606 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 607 if (snes->domainerror) { 608 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 609 PetscFunctionReturn(0); 610 } 611 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 612 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 613 if (neP->monitor) { 614 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr); 615 } 616 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */ 617 if (neP->monitor) { 618 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr); 619 } 620 goto theend1; 621 } 622 623 /* Fit points with cubic */ 624 count = 1; 625 while (PETSC_TRUE) { 626 if (lambda <= minlambda) { 627 if (neP->monitor) { 628 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); 629 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); 630 } 631 *flag = PETSC_FALSE; 632 break; 633 } 634 t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope; 635 t2 = .5*(gnormprev*gnormprev - fnorm*fnorm) - lambdaprev*initslope; 636 a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 637 b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); 638 d = b*b - 3*a*initslope; 639 if (d < 0.0) d = 0.0; 640 if (a == 0.0) { 641 lambdatemp = -initslope/(2.0*b); 642 } else { 643 lambdatemp = (-b + sqrt(d))/(3.0*a); 644 } 645 lambdaprev = lambda; 646 gnormprev = *gnorm; 647 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 648 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 649 else lambda = lambdatemp; 650 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 651 if (snes->nfuncs >= snes->max_funcs) { 652 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 653 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); 654 *flag = PETSC_FALSE; 655 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 656 break; 657 } 658 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 659 if (snes->domainerror) { 660 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 661 PetscFunctionReturn(0); 662 } 663 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 664 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 665 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* is reduction enough? */ 666 if (neP->monitor) { 667 ierr = PetscPrintf(comm," Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 668 } 669 break; 670 } else { 671 if (neP->monitor) { 672 ierr = PetscPrintf(comm," Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr); 673 } 674 } 675 count++; 676 } 677 theend1: 678 /* Optional user-defined check for line search step validity */ 679 if (neP->postcheckstep && *flag) { 680 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 681 if (changed_y) { 682 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 683 } 684 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 685 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 686 if (snes->domainerror) { 687 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 688 PetscFunctionReturn(0); 689 } 690 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 691 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 692 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 693 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 694 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 695 } 696 } 697 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 698 PetscFunctionReturn(0); 699 } 700 /* -------------------------------------------------------------------------- */ 701 #undef __FUNCT__ 702 #define __FUNCT__ "SNESLineSearchQuadratic" 703 /*@C 704 SNESLineSearchQuadratic - Performs a quadratic line search. 705 706 Collective on SNES and Vec 707 708 Input Parameters: 709 + snes - the SNES context 710 . lsctx - optional context for line search (not used here) 711 . x - current iterate 712 . f - residual evaluated at x 713 . y - search direction 714 . w - work vector 715 . fnorm - 2-norm of f 716 - xnorm - norm of x if known, otherwise 0 717 718 Output Parameters: 719 + g - residual evaluated at new iterate w 720 . w - new iterate (x + lambda*y) 721 . gnorm - 2-norm of g 722 . ynorm - 2-norm of search length 723 - flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure. 724 725 Options Database Keys: 726 + -snes_ls quadratic - Activates SNESLineSearchQuadratic() 727 . -snes_ls_alpha <alpha> - Sets alpha 728 . -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) 729 - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] ) 730 731 Notes: 732 Use SNESLineSearchSet() to set this routine within the SNESLS method. 733 734 Level: advanced 735 736 .keywords: SNES, nonlinear, quadratic, line search 737 738 .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms() 739 @*/ 740 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) 741 { 742 /* 743 Note that for line search purposes we work with with the related 744 minimization problem: 745 min z(x): R^n -> R, 746 where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2. 747 */ 748 PetscReal initslope,minlambda,lambda,lambdatemp,rellength; 749 #if defined(PETSC_USE_COMPLEX) 750 PetscScalar cinitslope; 751 #endif 752 PetscErrorCode ierr; 753 PetscInt count; 754 SNES_LS *neP = (SNES_LS*)snes->data; 755 PetscBool changed_w = PETSC_FALSE,changed_y = PETSC_FALSE; 756 757 PetscFunctionBegin; 758 ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 759 *flag = PETSC_TRUE; 760 761 ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr); 762 if (*ynorm == 0.0) { 763 if (neP->monitor) { 764 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr); 765 } 766 *gnorm = fnorm; 767 ierr = VecCopy(x,w);CHKERRQ(ierr); 768 ierr = VecCopy(f,g);CHKERRQ(ierr); 769 *flag = PETSC_FALSE; 770 goto theend2; 771 } 772 if (*ynorm > neP->maxstep) { /* Step too big, so scale back */ 773 ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr); 774 *ynorm = neP->maxstep; 775 } 776 ierr = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr); 777 minlambda = neP->minlambda/rellength; 778 ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr); 779 #if defined(PETSC_USE_COMPLEX) 780 ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr); 781 initslope = PetscRealPart(cinitslope); 782 #else 783 ierr = VecDot(f,w,&initslope);CHKERRQ(ierr); 784 #endif 785 if (initslope > 0.0) initslope = -initslope; 786 if (initslope == 0.0) initslope = -1.0; 787 ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr); 788 789 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 790 if (snes->nfuncs >= snes->max_funcs) { 791 ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); 792 *flag = PETSC_FALSE; 793 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 794 goto theend2; 795 } 796 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 797 if (snes->domainerror) { 798 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 799 PetscFunctionReturn(0); 800 } 801 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 802 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 803 if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */ 804 if (neP->monitor) { 805 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line Search: Using full step\n");CHKERRQ(ierr); 806 } 807 goto theend2; 808 } 809 810 /* Fit points with quadratic */ 811 lambda = 1.0; 812 count = 1; 813 while (PETSC_TRUE) { 814 if (lambda <= minlambda) { /* bad luck; use full step */ 815 if (neP->monitor) { 816 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr); 817 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); 818 } 819 ierr = VecCopy(x,w);CHKERRQ(ierr); 820 *flag = PETSC_FALSE; 821 break; 822 } 823 lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope); 824 if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; 825 if (lambdatemp <= .1*lambda) lambda = .1*lambda; 826 else lambda = lambdatemp; 827 828 ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr); 829 if (snes->nfuncs >= snes->max_funcs) { 830 ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); 831 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); 832 *flag = PETSC_FALSE; 833 snes->reason = SNES_DIVERGED_FUNCTION_COUNT; 834 break; 835 } 836 ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr); 837 if (snes->domainerror) { 838 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 839 PetscFunctionReturn(0); 840 } 841 ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr); 842 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 843 if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */ 844 if (neP->monitor) { 845 ierr = PetscViewerASCIIMonitorPrintf(neP->monitor," Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr); 846 } 847 break; 848 } 849 count++; 850 } 851 theend2: 852 /* Optional user-defined check for line search step validity */ 853 if (neP->postcheckstep) { 854 ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr); 855 if (changed_y) { 856 ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr); 857 } 858 if (changed_y || changed_w) { /* recompute the function if the step has changed */ 859 ierr = SNESComputeFunction(snes,w,g); 860 if (snes->domainerror) { 861 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 862 PetscFunctionReturn(0); 863 } 864 ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr); 865 ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr); 866 ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr); 867 ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr); 868 if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number"); 869 } 870 } 871 ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr); 872 PetscFunctionReturn(0); 873 } 874 875 /* -------------------------------------------------------------------------- */ 876 #undef __FUNCT__ 877 #define __FUNCT__ "SNESLineSearchSet" 878 /*@C 879 SNESLineSearchSet - Sets the line search routine to be used 880 by the method SNESLS. 881 882 Input Parameters: 883 + snes - nonlinear context obtained from SNESCreate() 884 . lsctx - optional user-defined context for use by line search 885 - func - pointer to int function 886 887 Logically Collective on SNES 888 889 Available Routines: 890 + SNESLineSearchCubic() - default line search 891 . SNESLineSearchQuadratic() - quadratic line search 892 . SNESLineSearchNo() - the full Newton step (actually not a line search) 893 - SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel) 894 895 Options Database Keys: 896 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 897 . -snes_ls_alpha <alpha> - Sets alpha 898 . -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) 899 - -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 900 901 Calling sequence of func: 902 .vb 903 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) 904 .ve 905 906 Input parameters for func: 907 + snes - nonlinear context 908 . lsctx - optional user-defined context for line search 909 . x - current iterate 910 . f - residual evaluated at x 911 . y - search direction 912 - fnorm - 2-norm of f 913 914 Output parameters for func: 915 + g - residual evaluated at new iterate y 916 . w - new iterate 917 . gnorm - 2-norm of g 918 . ynorm - 2-norm of search length 919 - flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure. 920 921 Level: advanced 922 923 .keywords: SNES, nonlinear, set, line search, routine 924 925 .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(), 926 SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck() 927 @*/ 928 PetscErrorCode SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *),void *lsctx) 929 { 930 PetscErrorCode ierr; 931 932 PetscFunctionBegin; 933 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); 934 PetscFunctionReturn(0); 935 } 936 937 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*/ 938 /* -------------------------------------------------------------------------- */ 939 EXTERN_C_BEGIN 940 #undef __FUNCT__ 941 #define __FUNCT__ "SNESLineSearchSet_LS" 942 PetscErrorCode SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx) 943 { 944 PetscFunctionBegin; 945 ((SNES_LS *)(snes->data))->LineSearch = func; 946 ((SNES_LS *)(snes->data))->lsP = lsctx; 947 PetscFunctionReturn(0); 948 } 949 EXTERN_C_END 950 /* -------------------------------------------------------------------------- */ 951 #undef __FUNCT__ 952 #define __FUNCT__ "SNESLineSearchSetPostCheck" 953 /*@C 954 SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed 955 by the line search routine in the Newton-based method SNESLS. 956 957 Input Parameters: 958 + snes - nonlinear context obtained from SNESCreate() 959 . func - pointer to function 960 - checkctx - optional user-defined context for use by step checking routine 961 962 Logically Collective on SNES 963 964 Calling sequence of func: 965 .vb 966 int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscBool *changed_y,PetscBool *changed_w) 967 .ve 968 where func returns an error code of 0 on success and a nonzero 969 on failure. 970 971 Input parameters for func: 972 + snes - nonlinear context 973 . checkctx - optional user-defined context for use by step checking routine 974 . x - previous iterate 975 . y - new search direction and length 976 - w - current candidate iterate 977 978 Output parameters for func: 979 + y - search direction (possibly changed) 980 . w - current iterate (possibly modified) 981 . changed_y - indicates search direction was changed by this routine 982 - changed_w - indicates current iterate was changed by this routine 983 984 Level: advanced 985 986 Notes: All line searches accept the new iterate computed by the line search checking routine. 987 988 Only one of changed_y and changed_w can be PETSC_TRUE 989 990 On input w = x + y 991 992 SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control 993 to the checking routine, and then (3) compute the corresponding nonlinear 994 function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}. 995 996 SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a 997 candidate nonlinear function f(u_{i+1}), (2) pass control to the checking 998 routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes 999 were made to the candidate iterate in the checking routine (as indicated 1000 by flag=PETSC_TRUE). The overhead of this extra function re-evaluation can be 1001 very costly, so use this feature with caution! 1002 1003 .keywords: SNES, nonlinear, set, line search check, step check, routine 1004 1005 .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck() 1006 @*/ 1007 PetscErrorCode SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void *checkctx) 1008 { 1009 PetscErrorCode ierr; 1010 1011 PetscFunctionBegin; 1012 ierr = PetscTryMethod(snes,"SNESLineSearchSetPostCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void*),(snes,func,checkctx));CHKERRQ(ierr); 1013 PetscFunctionReturn(0); 1014 } 1015 1016 #undef __FUNCT__ 1017 #define __FUNCT__ "SNESLineSearchSetPreCheck" 1018 /*@C 1019 SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve 1020 before the line search is called. 1021 1022 Input Parameters: 1023 + snes - nonlinear context obtained from SNESCreate() 1024 . func - pointer to function 1025 - checkctx - optional user-defined context for use by step checking routine 1026 1027 Logically Collective on SNES 1028 1029 Calling sequence of func: 1030 .vb 1031 int func (SNES snes, Vec x,Vec y,void *checkctx, PetscBool *changed_y) 1032 .ve 1033 where func returns an error code of 0 on success and a nonzero 1034 on failure. 1035 1036 Input parameters for func: 1037 + snes - nonlinear context 1038 . checkctx - optional user-defined context for use by step checking routine 1039 . x - previous iterate 1040 - y - new search direction and length 1041 1042 Output parameters for func: 1043 + y - search direction (possibly changed) 1044 - changed_y - indicates search direction was changed by this routine 1045 1046 Level: advanced 1047 1048 Notes: All line searches accept the new iterate computed by the line search checking routine. 1049 1050 .keywords: SNES, nonlinear, set, line search check, step check, routine 1051 1052 .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 1053 @*/ 1054 PetscErrorCode SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscBool *),void *checkctx) 1055 { 1056 PetscErrorCode ierr; 1057 1058 PetscFunctionBegin; 1059 ierr = PetscTryMethod(snes,"SNESLineSearchSetPreCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscBool *),void*),(snes,func,checkctx));CHKERRQ(ierr); 1060 PetscFunctionReturn(0); 1061 } 1062 1063 #undef __FUNCT__ 1064 #define __FUNCT__ "SNESLineSearchSetMonitor" 1065 /*@C 1066 SNESLineSearchSetMonitor - Prints information about the progress or lack of progress of the line search 1067 1068 Input Parameters: 1069 + snes - nonlinear context obtained from SNESCreate() 1070 - flg - PETSC_TRUE to monitor the line search 1071 1072 Logically Collective on SNES 1073 1074 Options Database: 1075 . -snes_ls_monitor 1076 1077 Level: intermediate 1078 1079 1080 .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate() 1081 @*/ 1082 PetscErrorCode SNESLineSearchSetMonitor(SNES snes,PetscBool flg) 1083 { 1084 PetscErrorCode ierr; 1085 1086 PetscFunctionBegin; 1087 ierr = PetscTryMethod(snes,"SNESLineSearchSetMonitor_C",(SNES,PetscBool),(snes,flg));CHKERRQ(ierr); 1088 PetscFunctionReturn(0); 1089 } 1090 1091 /* -------------------------------------------------------------------------- */ 1092 typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *); /* force argument to next function to not be extern C*/ 1093 typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscBool *); /* force argument to next function to not be extern C*/ 1094 EXTERN_C_BEGIN 1095 #undef __FUNCT__ 1096 #define __FUNCT__ "SNESLineSearchSetPostCheck_LS" 1097 PetscErrorCode SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx) 1098 { 1099 PetscFunctionBegin; 1100 ((SNES_LS *)(snes->data))->postcheckstep = func; 1101 ((SNES_LS *)(snes->data))->postcheck = checkctx; 1102 PetscFunctionReturn(0); 1103 } 1104 EXTERN_C_END 1105 1106 EXTERN_C_BEGIN 1107 #undef __FUNCT__ 1108 #define __FUNCT__ "SNESLineSearchSetPreCheck_LS" 1109 PetscErrorCode SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx) 1110 { 1111 PetscFunctionBegin; 1112 ((SNES_LS *)(snes->data))->precheckstep = func; 1113 ((SNES_LS *)(snes->data))->precheck = checkctx; 1114 PetscFunctionReturn(0); 1115 } 1116 EXTERN_C_END 1117 1118 EXTERN_C_BEGIN 1119 #undef __FUNCT__ 1120 #define __FUNCT__ "SNESLineSearchSetMonitor_LS" 1121 PetscErrorCode SNESLineSearchSetMonitor_LS(SNES snes,PetscBool flg) 1122 { 1123 SNES_LS *ls = (SNES_LS*)snes->data; 1124 PetscErrorCode ierr; 1125 1126 PetscFunctionBegin; 1127 if (flg && !ls->monitor) { 1128 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&ls->monitor);CHKERRQ(ierr); 1129 } else if (!flg && ls->monitor) { 1130 ierr = PetscViewerASCIIMonitorDestroy(ls->monitor);CHKERRQ(ierr); 1131 } 1132 PetscFunctionReturn(0); 1133 } 1134 EXTERN_C_END 1135 1136 /* 1137 SNESView_LS - Prints info from the SNESLS data structure. 1138 1139 Input Parameters: 1140 . SNES - the SNES context 1141 . viewer - visualization context 1142 1143 Application Interface Routine: SNESView() 1144 */ 1145 #undef __FUNCT__ 1146 #define __FUNCT__ "SNESView_LS" 1147 static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer) 1148 { 1149 SNES_LS *ls = (SNES_LS *)snes->data; 1150 const char *cstr; 1151 PetscErrorCode ierr; 1152 PetscBool iascii; 1153 1154 PetscFunctionBegin; 1155 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1156 if (iascii) { 1157 if (ls->LineSearch == SNESLineSearchNo) cstr = "SNESLineSearchNo"; 1158 else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic"; 1159 else if (ls->LineSearch == SNESLineSearchCubic) cstr = "SNESLineSearchCubic"; 1160 else cstr = "unknown"; 1161 ierr = PetscViewerASCIIPrintf(viewer," line search variant: %s\n",cstr);CHKERRQ(ierr); 1162 ierr = PetscViewerASCIIPrintf(viewer," alpha=%G, maxstep=%G, minlambda=%G\n",ls->alpha,ls->maxstep,ls->minlambda);CHKERRQ(ierr); 1163 } else { 1164 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name); 1165 } 1166 PetscFunctionReturn(0); 1167 } 1168 /* -------------------------------------------------------------------------- */ 1169 /* 1170 SNESSetFromOptions_LS - Sets various parameters for the SNESLS method. 1171 1172 Input Parameter: 1173 . snes - the SNES context 1174 1175 Application Interface Routine: SNESSetFromOptions() 1176 */ 1177 #undef __FUNCT__ 1178 #define __FUNCT__ "SNESSetFromOptions_LS" 1179 static PetscErrorCode SNESSetFromOptions_LS(SNES snes) 1180 { 1181 SNES_LS *ls = (SNES_LS *)snes->data; 1182 const char *lses[] = {"basic","basicnonorms","quadratic","cubic"}; 1183 PetscErrorCode ierr; 1184 PetscInt indx; 1185 PetscBool flg,set; 1186 1187 PetscFunctionBegin; 1188 ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr); 1189 ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr); 1190 ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr); 1191 ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",ls->minlambda,&ls->minlambda,0);CHKERRQ(ierr); 1192 ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",ls->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 1193 if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);} 1194 1195 ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr); 1196 if (flg) { 1197 switch (indx) { 1198 case 0: 1199 ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 1200 break; 1201 case 1: 1202 ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 1203 break; 1204 case 2: 1205 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr); 1206 break; 1207 case 3: 1208 ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr); 1209 break; 1210 } 1211 } 1212 ierr = PetscOptionsTail();CHKERRQ(ierr); 1213 PetscFunctionReturn(0); 1214 } 1215 /* -------------------------------------------------------------------------- */ 1216 /*MC 1217 SNESLS - Newton based nonlinear solver that uses a line search 1218 1219 Options Database: 1220 + -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search 1221 . -snes_ls_alpha <alpha> - Sets alpha 1222 . -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) 1223 . -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda / max_i ( y[i]/x[i] ) 1224 - -snes_ls_monitor - print information about progress of line searches 1225 1226 1227 Notes: This is the default nonlinear solver in SNES 1228 1229 Level: beginner 1230 1231 .seealso: SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(), 1232 SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(), 1233 SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams() 1234 1235 M*/ 1236 EXTERN_C_BEGIN 1237 #undef __FUNCT__ 1238 #define __FUNCT__ "SNESCreate_LS" 1239 PetscErrorCode SNESCreate_LS(SNES snes) 1240 { 1241 PetscErrorCode ierr; 1242 SNES_LS *neP; 1243 1244 PetscFunctionBegin; 1245 snes->ops->setup = SNESSetUp_LS; 1246 snes->ops->solve = SNESSolve_LS; 1247 snes->ops->destroy = SNESDestroy_LS; 1248 snes->ops->setfromoptions = SNESSetFromOptions_LS; 1249 snes->ops->view = SNESView_LS; 1250 1251 ierr = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr); 1252 snes->data = (void*)neP; 1253 neP->alpha = 1.e-4; 1254 neP->maxstep = 1.e8; 1255 neP->minlambda = 1.e-12; 1256 neP->LineSearch = SNESLineSearchCubic; 1257 neP->lsP = PETSC_NULL; 1258 neP->postcheckstep = PETSC_NULL; 1259 neP->postcheck = PETSC_NULL; 1260 neP->precheckstep = PETSC_NULL; 1261 neP->precheck = PETSC_NULL; 1262 1263 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_LS",SNESLineSearchSetMonitor_LS);CHKERRQ(ierr); 1264 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);CHKERRQ(ierr); 1265 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr); 1266 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr); 1267 1268 PetscFunctionReturn(0); 1269 } 1270 EXTERN_C_END 1271 1272 1273 1274