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