xref: /petsc/src/snes/linesearch/impls/cp/linesearchcp.c (revision 503c0ea9b45bcfbcebbb1ea5341243bbc69f0bea)
1 #include <petsc/private/linesearchimpl.h>
2 #include <petscsnes.h>
3 
4 static PetscErrorCode SNESLineSearchApply_CP(SNESLineSearch linesearch)
5 {
6   PetscBool      changed_y, changed_w;
7   Vec            X, Y, F, W;
8   SNES           snes;
9   PetscReal      xnorm, ynorm, gnorm, steptol, atol, rtol, ltol, maxstep;
10   PetscReal      lambda, lambda_old, lambda_update, delLambda;
11   PetscScalar    fty, fty_init, fty_old, fty_mid1, fty_mid2, s;
12   PetscInt       i, max_its;
13   PetscViewer    monitor;
14 
15   PetscFunctionBegin;
16   PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL));
17   PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm));
18   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
19   PetscCall(SNESLineSearchGetLambda(linesearch, &lambda));
20   PetscCall(SNESLineSearchGetTolerances(linesearch, &steptol, &maxstep, &rtol, &atol, &ltol, &max_its));
21   PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
22   PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
23 
24   /* precheck */
25   PetscCall(SNESLineSearchPreCheck(linesearch,X,Y,&changed_y));
26   lambda_old = 0.0;
27 
28   PetscCall(VecDot(F,Y,&fty_old));
29   if (PetscAbsScalar(fty_old) < atol * ynorm) {
30     if (monitor) {
31       PetscCall(PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel));
32       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));
33       PetscCall(PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel));
34     }
35     PetscCall(SNESSetConvergedReason(linesearch->snes,SNES_CONVERGED_FNORM_ABS));
36     PetscFunctionReturn(0);
37   }
38 
39   fty_init = fty_old;
40 
41   for (i = 0; i < max_its; i++) {
42     /* compute the norm at lambda */
43     PetscCall(VecCopy(X, W));
44     PetscCall(VecAXPY(W, -lambda, Y));
45     if (linesearch->ops->viproject) {
46       PetscCall((*linesearch->ops->viproject)(snes, W));
47     }
48     PetscCall((*linesearch->ops->snesfunc)(snes,W,F));
49     PetscCall(VecDot(F,Y,&fty));
50 
51     delLambda = lambda - lambda_old;
52 
53     /* check for convergence */
54     if (PetscAbsReal(delLambda) < steptol*lambda) break;
55     if (PetscAbsScalar(fty) / PetscAbsScalar(fty_init) < rtol) break;
56     if (PetscAbsScalar(fty) < atol * ynorm && i > 0) break;
57     if (monitor) {
58       PetscCall(PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel));
59       PetscCall(PetscViewerASCIIPrintf(monitor,"    Line search: lambdas = [%g, %g], ftys = [%g, %g]\n",(double)lambda, (double)lambda_old, (double)PetscRealPart(fty), (double)PetscRealPart(fty_old)));
60       PetscCall(PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel));
61     }
62 
63     /* compute the search direction */
64     if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) {
65       s = (fty - fty_old) / delLambda;
66     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
67       PetscCall(VecCopy(X, W));
68       PetscCall(VecAXPY(W, -0.5*(lambda + lambda_old), Y));
69       if (linesearch->ops->viproject) {
70         PetscCall((*linesearch->ops->viproject)(snes, W));
71       }
72       PetscCall((*linesearch->ops->snesfunc)(snes,W,F));
73       PetscCall(VecDot(F, Y, &fty_mid1));
74       s    = (3.*fty - 4.*fty_mid1 + fty_old) / delLambda;
75     } else {
76       PetscCall(VecCopy(X, W));
77       PetscCall(VecAXPY(W, -0.5*(lambda + lambda_old), Y));
78       if (linesearch->ops->viproject) {
79         PetscCall((*linesearch->ops->viproject)(snes, W));
80       }
81       PetscCall((*linesearch->ops->snesfunc)(snes,W,F));
82       PetscCall(VecDot(F, Y, &fty_mid1));
83       PetscCall(VecCopy(X, W));
84       PetscCall(VecAXPY(W, -(lambda + 0.5*(lambda - lambda_old)), Y));
85       if (linesearch->ops->viproject) {
86         PetscCall((*linesearch->ops->viproject)(snes, W));
87       }
88       PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
89       PetscCall(VecDot(F, Y, &fty_mid2));
90       s    = (2.*fty_mid2 + 3.*fty - 6.*fty_mid1 + fty_old) / (3.*delLambda);
91     }
92     /* if the solve is going in the wrong direction, fix it */
93     if (PetscRealPart(s) > 0.) s = -s;
94     if (s == 0.0) break;
95     lambda_update =  lambda - PetscRealPart(fty / s);
96 
97     /* switch directions if we stepped out of bounds */
98     if (lambda_update < steptol) lambda_update = lambda + PetscRealPart(fty / s);
99 
100     if (PetscIsInfOrNanReal(lambda_update)) break;
101     if (lambda_update > maxstep) break;
102 
103     /* compute the new state of the line search */
104     lambda_old = lambda;
105     lambda     = lambda_update;
106     fty_old    = fty;
107   }
108   /* construct the solution */
109   PetscCall(VecCopy(X, W));
110   PetscCall(VecAXPY(W, -lambda, Y));
111   if (linesearch->ops->viproject) {
112     PetscCall((*linesearch->ops->viproject)(snes, W));
113   }
114   /* postcheck */
115   PetscCall(SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w));
116   if (changed_y) {
117     PetscCall(VecAXPY(X, -lambda, Y));
118     if (linesearch->ops->viproject) {
119       PetscCall((*linesearch->ops->viproject)(snes, X));
120     }
121   } else {
122     PetscCall(VecCopy(W, X));
123   }
124   PetscCall((*linesearch->ops->snesfunc)(snes,X,F));
125 
126   PetscCall(SNESLineSearchComputeNorms(linesearch));
127   PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm));
128 
129   PetscCall(SNESLineSearchSetLambda(linesearch, lambda));
130 
131   if (monitor) {
132     PetscCall(PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel));
133     PetscCall(PetscViewerASCIIPrintf(monitor,"    Line search terminated: lambda = %g, fnorms = %g\n", (double)lambda, (double)gnorm));
134     PetscCall(PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel));
135   }
136   if (lambda <= steptol) {
137     PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT));
138   }
139   PetscFunctionReturn(0);
140 }
141 
142 /*MC
143    SNESLINESEARCHCP - Critical point line search. This line search assumes that there exists some
144    artificial G(x) for which the SNESFunction F(x) = grad G(x).  Therefore, this line search seeks
145    to find roots of dot(F, Y) via a secant method.
146 
147    Options Database Keys:
148 +  -snes_linesearch_minlambda <minlambda> - the minimum acceptable lambda
149 .  -snes_linesearch_maxstep <length> - the algorithm insures that a step length is never longer than this value
150 .  -snes_linesearch_damping <damping> - initial trial step length is scaled by this factor, default is 1.0
151 -  -snes_linesearch_max_it <max_it> - the maximum number of secant steps performed.
152 
153    Notes:
154    This method does NOT use the objective function if it is provided with SNESSetObjective().
155 
156    This method is the preferred line search for SNESQN and SNESNCG.
157 
158    Level: advanced
159 
160 .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
161 M*/
162 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch)
163 {
164   PetscFunctionBegin;
165   linesearch->ops->apply          = SNESLineSearchApply_CP;
166   linesearch->ops->destroy        = NULL;
167   linesearch->ops->setfromoptions = NULL;
168   linesearch->ops->reset          = NULL;
169   linesearch->ops->view           = NULL;
170   linesearch->ops->setup          = NULL;
171   linesearch->order               = SNES_LINESEARCH_ORDER_LINEAR;
172 
173   linesearch->max_its = 1;
174   PetscFunctionReturn(0);
175 }
176