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) { 25 PetscCall((*linesearch->ops->viproject)(snes, W)); 26 } 27 28 /* postcheck */ 29 PetscCall(SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w)); 30 if (changed_y) { 31 PetscCall(VecWAXPY(W,-lambda,Y,X)); 32 if (linesearch->ops->viproject) { 33 PetscCall((*linesearch->ops->viproject)(snes, W)); 34 } 35 } 36 if (linesearch->norms || snes->iter < snes->max_its-1) { 37 PetscCall((*linesearch->ops->snesfunc)(snes,W,F)); 38 PetscCall(SNESGetFunctionDomainError(snes, &domainerror)); 39 if (domainerror) { 40 PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_DOMAIN)); 41 PetscFunctionReturn(0); 42 } 43 } 44 45 if (linesearch->norms) { 46 if (!linesearch->ops->vinorm) PetscCall(VecNormBegin(F, NORM_2, &linesearch->fnorm)); 47 PetscCall(VecNormBegin(Y, NORM_2, &linesearch->ynorm)); 48 PetscCall(VecNormBegin(W, NORM_2, &linesearch->xnorm)); 49 if (!linesearch->ops->vinorm) PetscCall(VecNormEnd(F, NORM_2, &linesearch->fnorm)); 50 PetscCall(VecNormEnd(Y, NORM_2, &linesearch->ynorm)); 51 PetscCall(VecNormEnd(W, NORM_2, &linesearch->xnorm)); 52 53 if (linesearch->ops->vinorm) { 54 linesearch->fnorm = gnorm; 55 56 PetscCall((*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm)); 57 } else { 58 PetscCall(VecNorm(F,NORM_2,&linesearch->fnorm)); 59 } 60 } 61 62 /* copy the solution over */ 63 PetscCall(VecCopy(W, X)); 64 PetscFunctionReturn(0); 65 } 66 67 /*MC 68 SNESLINESEARCHBASIC - This line search implementation is not a line 69 search at all; it simply uses the full step. Thus, this routine is intended 70 for methods with well-scaled updates; i.e. Newton's method (SNESNEWTONLS), on 71 well-behaved problems. Also named as SNESLINESEARCHNONE 72 73 Options Database Keys: 74 + -snes_linesearch_damping <damping> - search vector is scaled by this amount, default is 1.0 75 - -snes_linesearch_norms <flag> - whether to compute norms or not, default is true (SNESLineSearchSetComputeNorms()) 76 77 Notes: 78 For methods with ill-scaled updates (SNESNRICHARDSON, SNESNCG), a small 79 damping parameter may yield satisfactory but slow convergence despite 80 the simplicity of the line search. 81 82 Level: advanced 83 84 .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNESLineSearchSetDamping()`, `SNESLineSearchSetComputeNorms()` 85 M*/ 86 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch) 87 { 88 PetscFunctionBegin; 89 linesearch->ops->apply = SNESLineSearchApply_Basic; 90 linesearch->ops->destroy = NULL; 91 linesearch->ops->setfromoptions = NULL; 92 linesearch->ops->reset = NULL; 93 linesearch->ops->view = NULL; 94 linesearch->ops->setup = NULL; 95 PetscFunctionReturn(0); 96 } 97