xref: /petsc/src/snes/linesearch/impls/cp/linesearchcp.c (revision ef46b1a67e276116c83b5d4ce8efc2932ea4fc0a)
1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h>
26a388c36SPeter Brune #include <petscsnes.h>
3bbd5d0b3SPeter Brune 
4f1c6b773SPeter Brune static PetscErrorCode SNESLineSearchApply_CP(SNESLineSearch linesearch)
5bbd5d0b3SPeter Brune {
6422a814eSBarry Smith   PetscBool      changed_y, changed_w;
76a388c36SPeter Brune   Vec            X, Y, F, W;
86a388c36SPeter Brune   SNES           snes;
9516fe3c3SPeter Brune   PetscReal      xnorm, ynorm, gnorm, steptol, atol, rtol, ltol, maxstep;
10bbd5d0b3SPeter Brune   PetscReal      lambda, lambda_old, lambda_update, delLambda;
11b8ac21a2SPeter Brune   PetscScalar    fty, fty_init, fty_old, fty_mid1, fty_mid2, s;
126a388c36SPeter Brune   PetscInt       i, max_its;
136a388c36SPeter Brune   PetscViewer    monitor;
14bbd5d0b3SPeter Brune 
15bbd5d0b3SPeter Brune   PetscFunctionBegin;
169566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL));
179566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm));
189566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
199566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetLambda(linesearch, &lambda));
209566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetTolerances(linesearch, &steptol, &maxstep, &rtol, &atol, &ltol, &max_its));
219566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
229566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
236a388c36SPeter Brune 
24bbd5d0b3SPeter Brune   /* precheck */
259566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchPreCheck(linesearch,X,Y,&changed_y));
26bbd5d0b3SPeter Brune   lambda_old = 0.0;
2762893cf2SPeter Brune 
289566063dSJacob Faibussowitsch   PetscCall(VecDot(F,Y,&fty_old));
2996b7d5c2SJed Brown   if (PetscAbsScalar(fty_old) < atol * ynorm) {
309e0a2ef9SBarry Smith     if (monitor) {
319566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel));
329566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(monitor,"    Line search terminated at initial point because dot(F,Y) = %g < atol*||y|| = %g\n",(double)PetscAbsScalar(fty_old), (double)atol*ynorm));
339566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel));
349e0a2ef9SBarry Smith     }
359566063dSJacob Faibussowitsch     PetscCall(SNESSetConvergedReason(linesearch->snes,SNES_CONVERGED_FNORM_ABS));
369e0a2ef9SBarry Smith     PetscFunctionReturn(0);
379e0a2ef9SBarry Smith   }
389e0a2ef9SBarry Smith 
3962893cf2SPeter Brune   fty_init = fty_old;
40c0b2db79SJed Brown 
41c0b2db79SJed Brown   for (i = 0; i < max_its; i++) {
42bbd5d0b3SPeter Brune     /* compute the norm at lambda */
43*ef46b1a6SStefano Zampini     PetscCall(VecWAXPY(W, -lambda, Y, X));
449bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
459566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->viproject)(snes, W));
469bd66eb0SPeter Brune     }
479566063dSJacob Faibussowitsch     PetscCall((*linesearch->ops->snesfunc)(snes,W,F));
489566063dSJacob Faibussowitsch     PetscCall(VecDot(F,Y,&fty));
49bbd5d0b3SPeter Brune 
50bbd5d0b3SPeter Brune     delLambda = lambda - lambda_old;
51b8ac21a2SPeter Brune 
52b8ac21a2SPeter Brune     /* check for convergence */
53b8ac21a2SPeter Brune     if (PetscAbsReal(delLambda) < steptol*lambda) break;
54b8ac21a2SPeter Brune     if (PetscAbsScalar(fty) / PetscAbsScalar(fty_init) < rtol) break;
5596b7d5c2SJed Brown     if (PetscAbsScalar(fty) < atol * ynorm && i > 0) break;
566a388c36SPeter Brune     if (monitor) {
579566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel));
589566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(monitor,"    Line search: lambdas = [%g, %g], ftys = [%g, %g]\n",(double)lambda, (double)lambda_old, (double)PetscRealPart(fty), (double)PetscRealPart(fty_old)));
599566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel));
60bbd5d0b3SPeter Brune     }
61bbd5d0b3SPeter Brune 
62bbd5d0b3SPeter Brune     /* compute the search direction */
63b000cd8dSPeter Brune     if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) {
64b8ac21a2SPeter Brune       s = (fty - fty_old) / delLambda;
65b000cd8dSPeter Brune     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
66*ef46b1a6SStefano Zampini       PetscCall(VecWAXPY(W, -0.5*(lambda + lambda_old), Y, X));
67b8ac21a2SPeter Brune       if (linesearch->ops->viproject) {
689566063dSJacob Faibussowitsch         PetscCall((*linesearch->ops->viproject)(snes, W));
69b8ac21a2SPeter Brune       }
709566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->snesfunc)(snes,W,F));
719566063dSJacob Faibussowitsch       PetscCall(VecDot(F, Y, &fty_mid1));
72b8ac21a2SPeter Brune       s    = (3.*fty - 4.*fty_mid1 + fty_old) / delLambda;
73b8ac21a2SPeter Brune     } else {
74*ef46b1a6SStefano Zampini       PetscCall(VecWAXPY(W, -0.5*(lambda + lambda_old), Y, X));
75b8ac21a2SPeter Brune       if (linesearch->ops->viproject) {
769566063dSJacob Faibussowitsch         PetscCall((*linesearch->ops->viproject)(snes, W));
77b8ac21a2SPeter Brune       }
789566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->snesfunc)(snes,W,F));
799566063dSJacob Faibussowitsch       PetscCall(VecDot(F, Y, &fty_mid1));
80*ef46b1a6SStefano Zampini       PetscCall(VecWAXPY(W, -(lambda + 0.5*(lambda - lambda_old)), Y, X));
81b8ac21a2SPeter Brune       if (linesearch->ops->viproject) {
829566063dSJacob Faibussowitsch         PetscCall((*linesearch->ops->viproject)(snes, W));
83b8ac21a2SPeter Brune       }
849566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
859566063dSJacob Faibussowitsch       PetscCall(VecDot(F, Y, &fty_mid2));
86b8ac21a2SPeter Brune       s    = (2.*fty_mid2 + 3.*fty - 6.*fty_mid1 + fty_old) / (3.*delLambda);
87b8ac21a2SPeter Brune     }
88b8ac21a2SPeter Brune     /* if the solve is going in the wrong direction, fix it */
89b8ac21a2SPeter Brune     if (PetscRealPart(s) > 0.) s = -s;
90f2f03862SBarry Smith     if (s == 0.0) break;
91b8ac21a2SPeter Brune     lambda_update =  lambda - PetscRealPart(fty / s);
92b8ac21a2SPeter Brune 
93b8ac21a2SPeter Brune     /* switch directions if we stepped out of bounds */
94f5af7f23SKarl Rupp     if (lambda_update < steptol) lambda_update = lambda + PetscRealPart(fty / s);
95b8ac21a2SPeter Brune 
96620fc8aaSSatish Balay     if (PetscIsInfOrNanReal(lambda_update)) break;
97f5af7f23SKarl Rupp     if (lambda_update > maxstep) break;
98bbd5d0b3SPeter Brune 
99bbd5d0b3SPeter Brune     /* compute the new state of the line search */
100bbd5d0b3SPeter Brune     lambda_old = lambda;
101bbd5d0b3SPeter Brune     lambda     = lambda_update;
102bbd5d0b3SPeter Brune     fty_old    = fty;
103bbd5d0b3SPeter Brune   }
104bbd5d0b3SPeter Brune   /* construct the solution */
105*ef46b1a6SStefano Zampini   PetscCall(VecWAXPY(W, -lambda, Y, X));
1069bd66eb0SPeter Brune   if (linesearch->ops->viproject) {
1079566063dSJacob Faibussowitsch     PetscCall((*linesearch->ops->viproject)(snes, W));
1089bd66eb0SPeter Brune   }
109bbd5d0b3SPeter Brune   /* postcheck */
1109566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w));
111bbd5d0b3SPeter Brune   if (changed_y) {
1129566063dSJacob Faibussowitsch     PetscCall(VecAXPY(X, -lambda, Y));
1139bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
1149566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->viproject)(snes, X));
1159bd66eb0SPeter Brune     }
116bbd5d0b3SPeter Brune   } else {
1179566063dSJacob Faibussowitsch     PetscCall(VecCopy(W, X));
118bbd5d0b3SPeter Brune   }
1199566063dSJacob Faibussowitsch   PetscCall((*linesearch->ops->snesfunc)(snes,X,F));
120bbd5d0b3SPeter Brune 
1219566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchComputeNorms(linesearch));
1229566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm));
123bbd5d0b3SPeter Brune 
1249566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetLambda(linesearch, lambda));
1256a388c36SPeter Brune 
1266a388c36SPeter Brune   if (monitor) {
1279566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel));
1289566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(monitor,"    Line search terminated: lambda = %g, fnorms = %g\n", (double)lambda, (double)gnorm));
1299566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel));
130bbd5d0b3SPeter Brune   }
1316a388c36SPeter Brune   if (lambda <= steptol) {
1329566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT));
133bbd5d0b3SPeter Brune   }
134bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
135bbd5d0b3SPeter Brune }
136bbd5d0b3SPeter Brune 
137954494b2SJed Brown /*MC
1381a4f838cSPeter Brune    SNESLINESEARCHCP - Critical point line search. This line search assumes that there exists some
139cd7522eaSPeter Brune    artificial G(x) for which the SNESFunction F(x) = grad G(x).  Therefore, this line search seeks
140db609ea7SPeter Brune    to find roots of dot(F, Y) via a secant method.
141cd7522eaSPeter Brune 
142cd7522eaSPeter Brune    Options Database Keys:
143eb23ba39SBarry Smith +  -snes_linesearch_minlambda <minlambda> - the minimum acceptable lambda
144eb23ba39SBarry Smith .  -snes_linesearch_maxstep <length> - the algorithm insures that a step length is never longer than this value
145eb23ba39SBarry Smith .  -snes_linesearch_damping <damping> - initial trial step length is scaled by this factor, default is 1.0
146eb23ba39SBarry Smith -  -snes_linesearch_max_it <max_it> - the maximum number of secant steps performed.
147cd7522eaSPeter Brune 
148cd7522eaSPeter Brune    Notes:
149eb23ba39SBarry Smith    This method does NOT use the objective function if it is provided with SNESSetObjective().
150eb23ba39SBarry Smith 
151cd7522eaSPeter Brune    This method is the preferred line search for SNESQN and SNESNCG.
152954494b2SJed Brown 
153954494b2SJed Brown    Level: advanced
154954494b2SJed Brown 
155f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
156954494b2SJed Brown M*/
1578cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch)
158bbd5d0b3SPeter Brune {
159bbd5d0b3SPeter Brune   PetscFunctionBegin;
160f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_CP;
1610298fd71SBarry Smith   linesearch->ops->destroy        = NULL;
1620298fd71SBarry Smith   linesearch->ops->setfromoptions = NULL;
1630298fd71SBarry Smith   linesearch->ops->reset          = NULL;
1640298fd71SBarry Smith   linesearch->ops->view           = NULL;
1650298fd71SBarry Smith   linesearch->ops->setup          = NULL;
166b000cd8dSPeter Brune   linesearch->order               = SNES_LINESEARCH_ORDER_LINEAR;
167a491bab8SPeter Brune 
168a491bab8SPeter Brune   linesearch->max_its = 1;
169bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
170bbd5d0b3SPeter Brune }
171