xref: /petsc/src/snes/linesearch/impls/cp/linesearchcp.c (revision 96b7d5c2ee3a706d3074402536cdfbeae64adfa4)
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;
7bbd5d0b3SPeter Brune   PetscErrorCode ierr;
86a388c36SPeter Brune   Vec            X, Y, F, W;
96a388c36SPeter Brune   SNES           snes;
10516fe3c3SPeter Brune   PetscReal      xnorm, ynorm, gnorm, steptol, atol, rtol, ltol, maxstep;
11bbd5d0b3SPeter Brune   PetscReal      lambda, lambda_old, lambda_update, delLambda;
12b8ac21a2SPeter Brune   PetscScalar    fty, fty_init, fty_old, fty_mid1, fty_mid2, s;
136a388c36SPeter Brune   PetscInt       i, max_its;
146a388c36SPeter Brune   PetscViewer    monitor;
15bbd5d0b3SPeter Brune 
16bbd5d0b3SPeter Brune   PetscFunctionBegin;
170298fd71SBarry Smith   ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL);CHKERRQ(ierr);
18f1c6b773SPeter Brune   ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr);
19f1c6b773SPeter Brune   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
20f1c6b773SPeter Brune   ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr);
21f1c6b773SPeter Brune   ierr = SNESLineSearchGetTolerances(linesearch, &steptol, &maxstep, &rtol, &atol, &ltol, &max_its);CHKERRQ(ierr);
22422a814eSBarry Smith   ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED);CHKERRQ(ierr);
23dcf2fd19SBarry Smith   ierr = SNESLineSearchGetDefaultMonitor(linesearch, &monitor);CHKERRQ(ierr);
246a388c36SPeter Brune 
25bbd5d0b3SPeter Brune   /* precheck */
267b1df9c1SPeter Brune   ierr       = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr);
27bbd5d0b3SPeter Brune   lambda_old = 0.0;
2862893cf2SPeter Brune 
29e0e0c2afSPeter Brune   ierr = VecDot(F,Y,&fty_old);CHKERRQ(ierr);
30*96b7d5c2SJed Brown   if (PetscAbsScalar(fty_old) < atol * ynorm) {
319e0a2ef9SBarry Smith     if (monitor) {
329e0a2ef9SBarry Smith       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
33*96b7d5c2SJed Brown       ierr = PetscViewerASCIIPrintf(monitor,"    Line search terminated at initial point because dot(F,Y) = %g < atol*||y|| = %g\n",(double)PetscAbsScalar(fty_old), (double)atol*ynorm);CHKERRQ(ierr);
349e0a2ef9SBarry Smith       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
359e0a2ef9SBarry Smith     }
36*96b7d5c2SJed Brown     ierr = SNESSetConvergedReason(linesearch->snes,SNES_CONVERGED_FNORM_ABS);CHKERRQ(ierr);
379e0a2ef9SBarry Smith     PetscFunctionReturn(0);
389e0a2ef9SBarry Smith   }
399e0a2ef9SBarry Smith 
4062893cf2SPeter Brune   fty_init = fty_old;
41c0b2db79SJed Brown 
42c0b2db79SJed Brown   for (i = 0; i < max_its; i++) {
43bbd5d0b3SPeter Brune     /* compute the norm at lambda */
44bbd5d0b3SPeter Brune     ierr = VecCopy(X, W);CHKERRQ(ierr);
45bbd5d0b3SPeter Brune     ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr);
469bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
479bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
489bd66eb0SPeter Brune     }
49fb8e56e0SPeter Brune     ierr = (*linesearch->ops->snesfunc)(snes,W,F);CHKERRQ(ierr);
50e0e0c2afSPeter Brune     ierr = VecDot(F,Y,&fty);CHKERRQ(ierr);
51bbd5d0b3SPeter Brune 
52bbd5d0b3SPeter Brune     delLambda = lambda - lambda_old;
53b8ac21a2SPeter Brune 
54b8ac21a2SPeter Brune     /* check for convergence */
55b8ac21a2SPeter Brune     if (PetscAbsReal(delLambda) < steptol*lambda) break;
56b8ac21a2SPeter Brune     if (PetscAbsScalar(fty) / PetscAbsScalar(fty_init) < rtol) break;
57*96b7d5c2SJed Brown     if (PetscAbsScalar(fty) < atol * ynorm && i > 0) break;
586a388c36SPeter Brune     if (monitor) {
596a388c36SPeter Brune       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
60c69d1a72SBarry Smith       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: lambdas = [%g, %g], ftys = [%g, %g]\n",(double)lambda, (double)lambda_old, (double)PetscRealPart(fty), (double)PetscRealPart(fty_old));CHKERRQ(ierr);
616a388c36SPeter Brune       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
62bbd5d0b3SPeter Brune     }
63bbd5d0b3SPeter Brune 
64bbd5d0b3SPeter Brune     /* compute the search direction */
65b000cd8dSPeter Brune     if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) {
66b8ac21a2SPeter Brune       s = (fty - fty_old) / delLambda;
67b000cd8dSPeter Brune     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
68b8ac21a2SPeter Brune       ierr = VecCopy(X, W);CHKERRQ(ierr);
69b8ac21a2SPeter Brune       ierr = VecAXPY(W, -0.5*(lambda + lambda_old), Y);CHKERRQ(ierr);
70b8ac21a2SPeter Brune       if (linesearch->ops->viproject) {
71b8ac21a2SPeter Brune         ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
72b8ac21a2SPeter Brune       }
73fb8e56e0SPeter Brune       ierr = (*linesearch->ops->snesfunc)(snes,W,F);CHKERRQ(ierr);
74b8ac21a2SPeter Brune       ierr = VecDot(F, Y, &fty_mid1);CHKERRQ(ierr);
75b8ac21a2SPeter Brune       s    = (3.*fty - 4.*fty_mid1 + fty_old) / delLambda;
76b8ac21a2SPeter Brune     } else {
77b8ac21a2SPeter Brune       ierr = VecCopy(X, W);CHKERRQ(ierr);
78b8ac21a2SPeter Brune       ierr = VecAXPY(W, -0.5*(lambda + lambda_old), Y);CHKERRQ(ierr);
79b8ac21a2SPeter Brune       if (linesearch->ops->viproject) {
80b8ac21a2SPeter Brune         ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
81b8ac21a2SPeter Brune       }
82fb8e56e0SPeter Brune       ierr = (*linesearch->ops->snesfunc)(snes,W,F);CHKERRQ(ierr);
83e0e0c2afSPeter Brune       ierr = VecDot(F, Y, &fty_mid1);CHKERRQ(ierr);
84b8ac21a2SPeter Brune       ierr = VecCopy(X, W);CHKERRQ(ierr);
85b8ac21a2SPeter Brune       ierr = VecAXPY(W, -(lambda + 0.5*(lambda - lambda_old)), Y);CHKERRQ(ierr);
86b8ac21a2SPeter Brune       if (linesearch->ops->viproject) {
87b8ac21a2SPeter Brune         ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
88b8ac21a2SPeter Brune       }
89fb8e56e0SPeter Brune       ierr = (*linesearch->ops->snesfunc)(snes, W, F);CHKERRQ(ierr);
90e0e0c2afSPeter Brune       ierr = VecDot(F, Y, &fty_mid2);CHKERRQ(ierr);
91b8ac21a2SPeter Brune       s    = (2.*fty_mid2 + 3.*fty - 6.*fty_mid1 + fty_old) / (3.*delLambda);
92b8ac21a2SPeter Brune     }
93b8ac21a2SPeter Brune     /* if the solve is going in the wrong direction, fix it */
94b8ac21a2SPeter Brune     if (PetscRealPart(s) > 0.) s = -s;
95f2f03862SBarry Smith     if (s == 0.0) break;
96b8ac21a2SPeter Brune     lambda_update =  lambda - PetscRealPart(fty / s);
97b8ac21a2SPeter Brune 
98b8ac21a2SPeter Brune     /* switch directions if we stepped out of bounds */
99f5af7f23SKarl Rupp     if (lambda_update < steptol) lambda_update = lambda + PetscRealPart(fty / s);
100b8ac21a2SPeter Brune 
101620fc8aaSSatish Balay     if (PetscIsInfOrNanReal(lambda_update)) break;
102f5af7f23SKarl Rupp     if (lambda_update > maxstep) break;
103bbd5d0b3SPeter Brune 
104bbd5d0b3SPeter Brune     /* compute the new state of the line search */
105bbd5d0b3SPeter Brune     lambda_old = lambda;
106bbd5d0b3SPeter Brune     lambda     = lambda_update;
107bbd5d0b3SPeter Brune     fty_old    = fty;
108bbd5d0b3SPeter Brune   }
109bbd5d0b3SPeter Brune   /* construct the solution */
110bbd5d0b3SPeter Brune   ierr = VecCopy(X, W);CHKERRQ(ierr);
111bbd5d0b3SPeter Brune   ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr);
1129bd66eb0SPeter Brune   if (linesearch->ops->viproject) {
1139bd66eb0SPeter Brune     ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
1149bd66eb0SPeter Brune   }
115bbd5d0b3SPeter Brune   /* postcheck */
1167b1df9c1SPeter Brune   ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr);
117bbd5d0b3SPeter Brune   if (changed_y) {
118bbd5d0b3SPeter Brune     ierr = VecAXPY(X, -lambda, Y);CHKERRQ(ierr);
1199bd66eb0SPeter Brune     if (linesearch->ops->viproject) {
1209bd66eb0SPeter Brune       ierr = (*linesearch->ops->viproject)(snes, X);CHKERRQ(ierr);
1219bd66eb0SPeter Brune     }
122bbd5d0b3SPeter Brune   } else {
123bbd5d0b3SPeter Brune     ierr = VecCopy(W, X);CHKERRQ(ierr);
124bbd5d0b3SPeter Brune   }
125fb8e56e0SPeter Brune   ierr = (*linesearch->ops->snesfunc)(snes,X,F);CHKERRQ(ierr);
126bbd5d0b3SPeter Brune 
127f1c6b773SPeter Brune   ierr = SNESLineSearchComputeNorms(linesearch);CHKERRQ(ierr);
128f1c6b773SPeter Brune   ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr);
129bbd5d0b3SPeter Brune 
130f1c6b773SPeter Brune   ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr);
1316a388c36SPeter Brune 
1326a388c36SPeter Brune   if (monitor) {
1336a388c36SPeter Brune     ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
134c69d1a72SBarry Smith     ierr = PetscViewerASCIIPrintf(monitor,"    Line search terminated: lambda = %g, fnorms = %g\n", (double)lambda, (double)gnorm);CHKERRQ(ierr);
1356a388c36SPeter Brune     ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
136bbd5d0b3SPeter Brune   }
1376a388c36SPeter Brune   if (lambda <= steptol) {
138e9b602ebSSatish Balay     ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);CHKERRQ(ierr);
139bbd5d0b3SPeter Brune   }
140bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
141bbd5d0b3SPeter Brune }
142bbd5d0b3SPeter Brune 
143954494b2SJed Brown /*MC
1441a4f838cSPeter Brune    SNESLINESEARCHCP - Critical point line search. This line search assumes that there exists some
145cd7522eaSPeter Brune    artificial G(x) for which the SNESFunction F(x) = grad G(x).  Therefore, this line search seeks
146db609ea7SPeter Brune    to find roots of dot(F, Y) via a secant method.
147cd7522eaSPeter Brune 
148cd7522eaSPeter Brune    Options Database Keys:
149eb23ba39SBarry Smith +  -snes_linesearch_minlambda <minlambda> - the minimum acceptable lambda
150eb23ba39SBarry Smith .  -snes_linesearch_maxstep <length> - the algorithm insures that a step length is never longer than this value
151eb23ba39SBarry Smith .  -snes_linesearch_damping <damping> - initial trial step length is scaled by this factor, default is 1.0
152eb23ba39SBarry Smith -  -snes_linesearch_max_it <max_it> - the maximum number of secant steps performed.
153cd7522eaSPeter Brune 
154cd7522eaSPeter Brune    Notes:
155eb23ba39SBarry Smith    This method does NOT use the objective function if it is provided with SNESSetObjective().
156eb23ba39SBarry Smith 
157cd7522eaSPeter Brune    This method is the preferred line search for SNESQN and SNESNCG.
158954494b2SJed Brown 
159954494b2SJed Brown    Level: advanced
160954494b2SJed Brown 
161f1c6b773SPeter Brune .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
162954494b2SJed Brown M*/
1638cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch)
164bbd5d0b3SPeter Brune {
165bbd5d0b3SPeter Brune   PetscFunctionBegin;
166f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_CP;
1670298fd71SBarry Smith   linesearch->ops->destroy        = NULL;
1680298fd71SBarry Smith   linesearch->ops->setfromoptions = NULL;
1690298fd71SBarry Smith   linesearch->ops->reset          = NULL;
1700298fd71SBarry Smith   linesearch->ops->view           = NULL;
1710298fd71SBarry Smith   linesearch->ops->setup          = NULL;
172b000cd8dSPeter Brune   linesearch->order               = SNES_LINESEARCH_ORDER_LINEAR;
173a491bab8SPeter Brune 
174a491bab8SPeter Brune   linesearch->max_its = 1;
175bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
176bbd5d0b3SPeter Brune }
177