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; 13*9566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL)); 14*9566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm)); 15*9566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetLambda(linesearch, &lambda)); 16*9566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetSNES(linesearch, &snes)); 17*9566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED)); 18bf7f4e0aSPeter Brune 19bf7f4e0aSPeter Brune /* precheck */ 20*9566063dSJacob Faibussowitsch PetscCall(SNESLineSearchPreCheck(linesearch,X,Y,&changed_y)); 21bf7f4e0aSPeter Brune 22bf7f4e0aSPeter Brune /* update */ 23*9566063dSJacob Faibussowitsch PetscCall(VecWAXPY(W,-lambda,Y,X)); 249bd66eb0SPeter Brune if (linesearch->ops->viproject) { 25*9566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->viproject)(snes, W)); 269bd66eb0SPeter Brune } 27bf7f4e0aSPeter Brune 28bf7f4e0aSPeter Brune /* postcheck */ 29*9566063dSJacob Faibussowitsch PetscCall(SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w)); 30bf7f4e0aSPeter Brune if (changed_y) { 31*9566063dSJacob Faibussowitsch PetscCall(VecWAXPY(W,-lambda,Y,X)); 329bd66eb0SPeter Brune if (linesearch->ops->viproject) { 33*9566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->viproject)(snes, W)); 349bd66eb0SPeter Brune } 35bf7f4e0aSPeter Brune } 365fc09b23SPeter Brune if (linesearch->norms || snes->iter < snes->max_its-1) { 37*9566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->snesfunc)(snes,W,F)); 38*9566063dSJacob Faibussowitsch PetscCall(SNESGetFunctionDomainError(snes, &domainerror)); 396a388c36SPeter Brune if (domainerror) { 40*9566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_DOMAIN)); 41bf7f4e0aSPeter Brune PetscFunctionReturn(0); 42bf7f4e0aSPeter Brune } 435fc09b23SPeter Brune } 446a388c36SPeter Brune 455fc09b23SPeter Brune if (linesearch->norms) { 46*9566063dSJacob Faibussowitsch if (!linesearch->ops->vinorm) PetscCall(VecNormBegin(F, NORM_2, &linesearch->fnorm)); 47*9566063dSJacob Faibussowitsch PetscCall(VecNormBegin(Y, NORM_2, &linesearch->ynorm)); 48*9566063dSJacob Faibussowitsch PetscCall(VecNormBegin(W, NORM_2, &linesearch->xnorm)); 49*9566063dSJacob Faibussowitsch if (!linesearch->ops->vinorm) PetscCall(VecNormEnd(F, NORM_2, &linesearch->fnorm)); 50*9566063dSJacob Faibussowitsch PetscCall(VecNormEnd(Y, NORM_2, &linesearch->ynorm)); 51*9566063dSJacob Faibussowitsch PetscCall(VecNormEnd(W, NORM_2, &linesearch->xnorm)); 525fc09b23SPeter Brune 539bd66eb0SPeter Brune if (linesearch->ops->vinorm) { 549bd66eb0SPeter Brune linesearch->fnorm = gnorm; 55f5af7f23SKarl Rupp 56*9566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->vinorm)(snes, F, W, &linesearch->fnorm)); 579bd66eb0SPeter Brune } else { 58*9566063dSJacob Faibussowitsch PetscCall(VecNorm(F,NORM_2,&linesearch->fnorm)); 599bd66eb0SPeter Brune } 605fc09b23SPeter Brune } 616a388c36SPeter Brune 62bf7f4e0aSPeter Brune /* copy the solution over */ 63*9566063dSJacob Faibussowitsch PetscCall(VecCopy(W, X)); 64bf7f4e0aSPeter Brune PetscFunctionReturn(0); 65bf7f4e0aSPeter Brune } 66bf7f4e0aSPeter Brune 67954494b2SJed Brown /*MC 681a4f838cSPeter Brune SNESLINESEARCHBASIC - This line search implementation is not a line 69cd7522eaSPeter Brune search at all; it simply uses the full step. Thus, this routine is intended 7004d7464bSBarry Smith for methods with well-scaled updates; i.e. Newton's method (SNESNEWTONLS), on 71cd7522eaSPeter Brune well-behaved problems. 72cd7522eaSPeter Brune 73cd7522eaSPeter Brune Options Database Keys: 747c720adbSBarry Smith + -snes_linesearch_damping <damping> - search vector is scaled by this amount, default is 1.0 751fe24845SBarry Smith - -snes_linesearch_norms <flag> - whether to compute norms or not, default is true (SNESLineSearchSetComputeNorms()) 76cd7522eaSPeter Brune 77cd7522eaSPeter Brune Notes: 78cd7522eaSPeter Brune For methods with ill-scaled updates (SNESNRICHARDSON, SNESNCG), a small 79cd7522eaSPeter Brune damping parameter may yield satisfactory but slow convergence despite 80cd7522eaSPeter Brune the simplicity of the line search. 81954494b2SJed Brown 82954494b2SJed Brown Level: advanced 83954494b2SJed Brown 841fe24845SBarry Smith .seealso: SNESLineSearchCreate(), SNESLineSearchSetType(), SNESLineSearchSetDamping(), SNESLineSearchSetComputeNorms() 85954494b2SJed Brown M*/ 868cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Basic(SNESLineSearch linesearch) 87bf7f4e0aSPeter Brune { 88bf7f4e0aSPeter Brune PetscFunctionBegin; 89f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_Basic; 900298fd71SBarry Smith linesearch->ops->destroy = NULL; 910298fd71SBarry Smith linesearch->ops->setfromoptions = NULL; 920298fd71SBarry Smith linesearch->ops->reset = NULL; 930298fd71SBarry Smith linesearch->ops->view = NULL; 940298fd71SBarry Smith linesearch->ops->setup = NULL; 95bf7f4e0aSPeter Brune PetscFunctionReturn(0); 96bf7f4e0aSPeter Brune } 97