1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> 2af0996ceSBarry Smith #include <petsc/private/snesimpl.h> 3bf7f4e0aSPeter Brune 4f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchApply_Basic(SNESLineSearch linesearch) 5bf7f4e0aSPeter Brune { 6bf7f4e0aSPeter Brune PetscBool changed_y, changed_w; 76a388c36SPeter Brune Vec X, F, Y, W; 86a388c36SPeter Brune SNES snes; 96a388c36SPeter Brune PetscReal gnorm, xnorm, ynorm, lambda; 106a388c36SPeter Brune PetscBool domainerror; 11bf7f4e0aSPeter Brune 12bf7f4e0aSPeter Brune PetscFunctionBegin; 139566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL)); 149566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm)); 159566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetLambda(linesearch, &lambda)); 169566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 179566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED)); 18bf7f4e0aSPeter Brune 19bf7f4e0aSPeter Brune /* precheck */ 209566063dSJacob Faibussowitsch PetscCall(SNESLineSearchPreCheck(linesearch,X,Y,&changed_y)); 21bf7f4e0aSPeter Brune 22bf7f4e0aSPeter Brune /* update */ 239566063dSJacob Faibussowitsch PetscCall(VecWAXPY(W,-lambda,Y,X)); 24*1baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 25bf7f4e0aSPeter Brune 26bf7f4e0aSPeter Brune /* postcheck */ 279566063dSJacob Faibussowitsch PetscCall(SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w)); 28bf7f4e0aSPeter Brune if (changed_y) { 299566063dSJacob Faibussowitsch PetscCall(VecWAXPY(W,-lambda,Y,X)); 30*1baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W)); 31bf7f4e0aSPeter Brune } 325fc09b23SPeter Brune if (linesearch->norms || snes->iter < snes->max_its-1) { 339566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->snesfunc)(snes,W,F)); 349566063dSJacob Faibussowitsch PetscCall(SNESGetFunctionDomainError(snes, &domainerror)); 356a388c36SPeter Brune if (domainerror) { 369566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_DOMAIN)); 37bf7f4e0aSPeter Brune PetscFunctionReturn(0); 38bf7f4e0aSPeter Brune } 395fc09b23SPeter Brune } 406a388c36SPeter Brune 415fc09b23SPeter Brune if (linesearch->norms) { 429566063dSJacob Faibussowitsch if (!linesearch->ops->vinorm) PetscCall(VecNormBegin(F, NORM_2, &linesearch->fnorm)); 439566063dSJacob Faibussowitsch PetscCall(VecNormBegin(Y, NORM_2, &linesearch->ynorm)); 449566063dSJacob Faibussowitsch PetscCall(VecNormBegin(W, NORM_2, &linesearch->xnorm)); 459566063dSJacob Faibussowitsch if (!linesearch->ops->vinorm) PetscCall(VecNormEnd(F, NORM_2, &linesearch->fnorm)); 469566063dSJacob Faibussowitsch PetscCall(VecNormEnd(Y, NORM_2, &linesearch->ynorm)); 479566063dSJacob Faibussowitsch PetscCall(VecNormEnd(W, NORM_2, &linesearch->xnorm)); 485fc09b23SPeter Brune 499bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 509bd66eb0SPeter Brune linesearch->fnorm = gnorm; 51f5af7f23SKarl Rupp 529566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm)); 539bd66eb0SPeter Brune } else { 549566063dSJacob Faibussowitsch PetscCall(VecNorm(F,NORM_2,&linesearch->fnorm)); 559bd66eb0SPeter Brune } 565fc09b23SPeter Brune } 576a388c36SPeter Brune 58bf7f4e0aSPeter Brune /* copy the solution over */ 599566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X)); 60bf7f4e0aSPeter Brune PetscFunctionReturn(0); 61bf7f4e0aSPeter Brune } 62bf7f4e0aSPeter Brune 63954494b2SJed Brown /*MC 641a4f838cSPeter Brune SNESLINESEARCHBASIC - This line search implementation is not a line 65cd7522eaSPeter Brune search at all; it simply uses the full step. Thus, this routine is intended 6604d7464bSBarry Smith for methods with well-scaled updates; i.e. Newton's method (SNESNEWTONLS), on 670b00b554SBarry Smith well-behaved problems. Also named as SNESLINESEARCHNONE 68cd7522eaSPeter Brune 69cd7522eaSPeter Brune Options Database Keys: 707c720adbSBarry Smith + -snes_linesearch_damping <damping> - search vector is scaled by this amount, default is 1.0 711fe24845SBarry Smith - -snes_linesearch_norms <flag> - whether to compute norms or not, default is true (SNESLineSearchSetComputeNorms()) 72cd7522eaSPeter Brune 73cd7522eaSPeter Brune Notes: 74cd7522eaSPeter Brune For methods with ill-scaled updates (SNESNRICHARDSON, SNESNCG), a small 75cd7522eaSPeter Brune damping parameter may yield satisfactory but slow convergence despite 76cd7522eaSPeter Brune the simplicity of the line search. 77954494b2SJed Brown 78954494b2SJed Brown Level: advanced 79954494b2SJed Brown 80db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNESLineSearchSetDamping()`, `SNESLineSearchSetComputeNorms()` 81954494b2SJed Brown M*/ 828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch) 83bf7f4e0aSPeter Brune { 84bf7f4e0aSPeter Brune PetscFunctionBegin; 85f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_Basic; 860298fd71SBarry Smith linesearch->ops->destroy = NULL; 870298fd71SBarry Smith linesearch->ops->setfromoptions = NULL; 880298fd71SBarry Smith linesearch->ops->reset = NULL; 890298fd71SBarry Smith linesearch->ops->view = NULL; 900298fd71SBarry Smith linesearch->ops->setup = NULL; 91bf7f4e0aSPeter Brune PetscFunctionReturn(0); 92bf7f4e0aSPeter Brune } 93