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