xref: /petsc/src/snes/linesearch/impls/cp/linesearchcp.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h>
26a388c36SPeter Brune #include <petscsnes.h>
3bbd5d0b3SPeter Brune 
49371c9d4SSatish Balay static PetscErrorCode SNESLineSearchApply_CP(SNESLineSearch linesearch) {
5422a814eSBarry Smith   PetscBool   changed_y, changed_w;
66a388c36SPeter Brune   Vec         X, Y, F, W;
76a388c36SPeter Brune   SNES        snes;
8516fe3c3SPeter Brune   PetscReal   xnorm, ynorm, gnorm, steptol, atol, rtol, ltol, maxstep;
9bbd5d0b3SPeter Brune   PetscReal   lambda, lambda_old, lambda_update, delLambda;
10b8ac21a2SPeter Brune   PetscScalar fty, fty_init, fty_old, fty_mid1, fty_mid2, s;
116a388c36SPeter Brune   PetscInt    i, max_its;
126a388c36SPeter Brune   PetscViewer monitor;
13bbd5d0b3SPeter Brune 
14bbd5d0b3SPeter Brune   PetscFunctionBegin;
159566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL));
169566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm));
179566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
189566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetLambda(linesearch, &lambda));
199566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetTolerances(linesearch, &steptol, &maxstep, &rtol, &atol, &ltol, &max_its));
209566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
219566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
226a388c36SPeter Brune 
23bbd5d0b3SPeter Brune   /* precheck */
249566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y));
25bbd5d0b3SPeter Brune   lambda_old = 0.0;
2662893cf2SPeter Brune 
279566063dSJacob Faibussowitsch   PetscCall(VecDot(F, Y, &fty_old));
2896b7d5c2SJed Brown   if (PetscAbsScalar(fty_old) < atol * ynorm) {
299e0a2ef9SBarry Smith     if (monitor) {
309566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
3163a3b9bcSJacob 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)));
329566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
339e0a2ef9SBarry Smith     }
349566063dSJacob Faibussowitsch     PetscCall(SNESSetConvergedReason(linesearch->snes, SNES_CONVERGED_FNORM_ABS));
359e0a2ef9SBarry Smith     PetscFunctionReturn(0);
369e0a2ef9SBarry Smith   }
379e0a2ef9SBarry Smith 
3862893cf2SPeter Brune   fty_init = fty_old;
39c0b2db79SJed Brown 
40c0b2db79SJed Brown   for (i = 0; i < max_its; i++) {
41bbd5d0b3SPeter Brune     /* compute the norm at lambda */
42ef46b1a6SStefano Zampini     PetscCall(VecWAXPY(W, -lambda, Y, X));
431baa6e33SBarry Smith     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
449566063dSJacob Faibussowitsch     PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
459566063dSJacob Faibussowitsch     PetscCall(VecDot(F, Y, &fty));
46bbd5d0b3SPeter Brune 
47bbd5d0b3SPeter Brune     delLambda = lambda - lambda_old;
48b8ac21a2SPeter Brune 
49b8ac21a2SPeter Brune     /* check for convergence */
50b8ac21a2SPeter Brune     if (PetscAbsReal(delLambda) < steptol * lambda) break;
51b8ac21a2SPeter Brune     if (PetscAbsScalar(fty) / PetscAbsScalar(fty_init) < rtol) break;
5296b7d5c2SJed Brown     if (PetscAbsScalar(fty) < atol * ynorm && i > 0) break;
536a388c36SPeter Brune     if (monitor) {
549566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
559566063dSJacob 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)));
569566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
57bbd5d0b3SPeter Brune     }
58bbd5d0b3SPeter Brune 
59bbd5d0b3SPeter Brune     /* compute the search direction */
60b000cd8dSPeter Brune     if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) {
61b8ac21a2SPeter Brune       s = (fty - fty_old) / delLambda;
62b000cd8dSPeter Brune     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
63ef46b1a6SStefano Zampini       PetscCall(VecWAXPY(W, -0.5 * (lambda + lambda_old), Y, X));
641baa6e33SBarry Smith       if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
659566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
669566063dSJacob Faibussowitsch       PetscCall(VecDot(F, Y, &fty_mid1));
67b8ac21a2SPeter Brune       s = (3. * fty - 4. * fty_mid1 + fty_old) / delLambda;
68b8ac21a2SPeter Brune     } else {
69ef46b1a6SStefano Zampini       PetscCall(VecWAXPY(W, -0.5 * (lambda + lambda_old), Y, X));
701baa6e33SBarry Smith       if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
719566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
729566063dSJacob Faibussowitsch       PetscCall(VecDot(F, Y, &fty_mid1));
73ef46b1a6SStefano Zampini       PetscCall(VecWAXPY(W, -(lambda + 0.5 * (lambda - lambda_old)), Y, X));
741baa6e33SBarry Smith       if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
759566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
769566063dSJacob Faibussowitsch       PetscCall(VecDot(F, Y, &fty_mid2));
77b8ac21a2SPeter Brune       s = (2. * fty_mid2 + 3. * fty - 6. * fty_mid1 + fty_old) / (3. * delLambda);
78b8ac21a2SPeter Brune     }
79b8ac21a2SPeter Brune     /* if the solve is going in the wrong direction, fix it */
80b8ac21a2SPeter Brune     if (PetscRealPart(s) > 0.) s = -s;
81f2f03862SBarry Smith     if (s == 0.0) break;
82b8ac21a2SPeter Brune     lambda_update = lambda - PetscRealPart(fty / s);
83b8ac21a2SPeter Brune 
84b8ac21a2SPeter Brune     /* switch directions if we stepped out of bounds */
85f5af7f23SKarl Rupp     if (lambda_update < steptol) lambda_update = lambda + PetscRealPart(fty / s);
86b8ac21a2SPeter Brune 
87620fc8aaSSatish Balay     if (PetscIsInfOrNanReal(lambda_update)) break;
88f5af7f23SKarl Rupp     if (lambda_update > maxstep) break;
89bbd5d0b3SPeter Brune 
90bbd5d0b3SPeter Brune     /* compute the new state of the line search */
91bbd5d0b3SPeter Brune     lambda_old = lambda;
92bbd5d0b3SPeter Brune     lambda     = lambda_update;
93bbd5d0b3SPeter Brune     fty_old    = fty;
94bbd5d0b3SPeter Brune   }
95bbd5d0b3SPeter Brune   /* construct the solution */
96ef46b1a6SStefano Zampini   PetscCall(VecWAXPY(W, -lambda, Y, X));
971baa6e33SBarry Smith   if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
98bbd5d0b3SPeter Brune   /* postcheck */
999566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w));
100bbd5d0b3SPeter Brune   if (changed_y) {
1019566063dSJacob Faibussowitsch     PetscCall(VecAXPY(X, -lambda, Y));
1021baa6e33SBarry Smith     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, X));
103bbd5d0b3SPeter Brune   } else {
1049566063dSJacob Faibussowitsch     PetscCall(VecCopy(W, X));
105bbd5d0b3SPeter Brune   }
1069566063dSJacob Faibussowitsch   PetscCall((*linesearch->ops->snesfunc)(snes, X, F));
107bbd5d0b3SPeter Brune 
1089566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchComputeNorms(linesearch));
1099566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm));
110bbd5d0b3SPeter Brune 
1119566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetLambda(linesearch, lambda));
1126a388c36SPeter Brune 
1136a388c36SPeter Brune   if (monitor) {
1149566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
1159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search terminated: lambda = %g, fnorms = %g\n", (double)lambda, (double)gnorm));
1169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
117bbd5d0b3SPeter Brune   }
118*48a46eb9SPierre Jolivet   if (lambda <= steptol) PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT));
119bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
120bbd5d0b3SPeter Brune }
121bbd5d0b3SPeter Brune 
122954494b2SJed Brown /*MC
1231a4f838cSPeter Brune    SNESLINESEARCHCP - Critical point line search. This line search assumes that there exists some
124cd7522eaSPeter Brune    artificial G(x) for which the SNESFunction F(x) = grad G(x).  Therefore, this line search seeks
125db609ea7SPeter Brune    to find roots of dot(F, Y) via a secant method.
126cd7522eaSPeter Brune 
127cd7522eaSPeter Brune    Options Database Keys:
128eb23ba39SBarry Smith +  -snes_linesearch_minlambda <minlambda> - the minimum acceptable lambda
129eb23ba39SBarry Smith .  -snes_linesearch_maxstep <length> - the algorithm insures that a step length is never longer than this value
130eb23ba39SBarry Smith .  -snes_linesearch_damping <damping> - initial trial step length is scaled by this factor, default is 1.0
131eb23ba39SBarry Smith -  -snes_linesearch_max_it <max_it> - the maximum number of secant steps performed.
132cd7522eaSPeter Brune 
133cd7522eaSPeter Brune    Notes:
134eb23ba39SBarry Smith    This method does NOT use the objective function if it is provided with SNESSetObjective().
135eb23ba39SBarry Smith 
136cd7522eaSPeter Brune    This method is the preferred line search for SNESQN and SNESNCG.
137954494b2SJed Brown 
138954494b2SJed Brown    Level: advanced
139954494b2SJed Brown 
140db781477SPatrick Sanan .seealso: `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
141954494b2SJed Brown M*/
1429371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch) {
143bbd5d0b3SPeter Brune   PetscFunctionBegin;
144f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_CP;
1450298fd71SBarry Smith   linesearch->ops->destroy        = NULL;
1460298fd71SBarry Smith   linesearch->ops->setfromoptions = NULL;
1470298fd71SBarry Smith   linesearch->ops->reset          = NULL;
1480298fd71SBarry Smith   linesearch->ops->view           = NULL;
1490298fd71SBarry Smith   linesearch->ops->setup          = NULL;
150b000cd8dSPeter Brune   linesearch->order               = SNES_LINESEARCH_ORDER_LINEAR;
151a491bab8SPeter Brune 
152a491bab8SPeter Brune   linesearch->max_its = 1;
153bbd5d0b3SPeter Brune   PetscFunctionReturn(0);
154bbd5d0b3SPeter Brune }
155