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