1 #include <petsc/private/linesearchimpl.h> 2 #include <petsc/private/snesimpl.h> 3 4 static PetscErrorCode SNESLineSearchApply_Basic(SNESLineSearch linesearch) 5 { 6 PetscBool changed_y, changed_w; 7 Vec X, F, Y, W; 8 SNES snes; 9 PetscReal gnorm, xnorm, ynorm, lambda; 10 PetscBool domainerror; 11 12 PetscFunctionBegin; 13 PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL)); 14 PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm)); 15 PetscCall(SNESLineSearchGetLambda(linesearch, &lambda)); 16 PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 17 PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED)); 18 19 /* precheck */ 20 PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y)); 21 22 /* update */ 23 PetscCall(VecWAXPY(W, -lambda, Y, X)); 24 if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 25 26 /* postcheck */ 27 PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w)); 28 if (changed_y) { 29 if (!changed_w) PetscCall(VecWAXPY(W, -lambda, Y, X)); 30 if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 31 } 32 if (linesearch->norms || snes->iter < snes->max_its - 1) { 33 PetscCall((*linesearch->ops->snesfunc)(snes, W, F)); 34 PetscCall(SNESGetFunctionDomainError(snes, &domainerror)); 35 if (domainerror) { 36 PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_DOMAIN)); 37 PetscFunctionReturn(PETSC_SUCCESS); 38 } 39 } 40 if (linesearch->norms) { 41 if (!linesearch->ops->vinorm) PetscCall(VecNormBegin(F, NORM_2, &linesearch->fnorm)); 42 PetscCall(VecNormBegin(Y, NORM_2, &linesearch->ynorm)); 43 PetscCall(VecNormBegin(W, NORM_2, &linesearch->xnorm)); 44 if (!linesearch->ops->vinorm) PetscCall(VecNormEnd(F, NORM_2, &linesearch->fnorm)); 45 PetscCall(VecNormEnd(Y, NORM_2, &linesearch->ynorm)); 46 PetscCall(VecNormEnd(W, NORM_2, &linesearch->xnorm)); 47 48 if (linesearch->ops->vinorm) { 49 linesearch->fnorm = gnorm; 50 51 PetscCall((*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm)); 52 } else { 53 PetscCall(VecNorm(F, NORM_2, &linesearch->fnorm)); 54 } 55 } 56 57 /* copy the solution over */ 58 PetscCall(VecCopy(W, X)); 59 PetscFunctionReturn(PETSC_SUCCESS); 60 } 61 62 /*MC 63 SNESLINESEARCHBASIC - This line search implementation is not a line 64 search at all; it simply uses the full step $x_{k+1} = x_k - \lambda Y_k$ with $\lambda=1$. 65 Alternatively, $\lambda$ can be configured to be a constant damping factor by setting `snes_linesearch_damping`. 66 Thus, this routine is intended for methods with well-scaled updates; i.e. Newton's method (`SNESNEWTONLS`), on 67 well-behaved problems. Also named `SNESLINESEARCHNONE`. 68 69 Options Database Keys: 70 + -snes_linesearch_damping <1.0> - step length is scaled by this factor 71 - -snes_linesearch_norms <true> - whether to compute norms or not (`SNESLineSearchSetComputeNorms()`) 72 73 Note: 74 For methods with ill-scaled updates (`SNESNRICHARDSON`, `SNESNCG`), a small 75 damping parameter may yield satisfactory, but slow convergence, despite 76 the lack of the line search. 77 78 Level: advanced 79 80 .seealso: [](ch_snes), `SNES`, `SNESLineSearch`, `SNESLineSearchType`, `SNESGetLineSearch()`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNESLineSearchSetDamping()`, `SNESLineSearchSetComputeNorms()` 81 M*/ 82 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch) 83 { 84 PetscFunctionBegin; 85 linesearch->ops->apply = SNESLineSearchApply_Basic; 86 linesearch->ops->destroy = NULL; 87 linesearch->ops->setfromoptions = NULL; 88 linesearch->ops->reset = NULL; 89 linesearch->ops->view = NULL; 90 linesearch->ops->setup = NULL; 91 PetscFunctionReturn(PETSC_SUCCESS); 92 } 93