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