xref: /petsc/src/snes/linesearch/impls/cp/linesearchcp.c (revision 82d48c8790e4c36748aef416be2d8358100e934d)
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   ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL);CHKERRQ(ierr);
22   ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr);
23   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
24   ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr);
25   ierr = SNESLineSearchGetTolerances(linesearch, &steptol, &maxstep, &rtol, &atol, &ltol, &max_its);CHKERRQ(ierr);
26   ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr);
27   ierr = SNESLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr);
28 
29   /* precheck */
30   ierr       = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr);
31   lambda_old = 0.0;
32 
33 
34   for (i = 0; i < max_its; i++) {
35     if (i == 0) {
36       ierr = VecDot(F,Y,&fty_old);CHKERRQ(ierr);
37       fty_init = fty_old;
38     }
39     /* compute the norm at lambda */
40     ierr = VecCopy(X, W);CHKERRQ(ierr);
41     ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr);
42     if (linesearch->ops->viproject) {
43       ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
44     }
45     ierr = (*linesearch->ops->snesfunc)(snes,W,F);CHKERRQ(ierr);
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 = (*linesearch->ops->snesfunc)(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 = (*linesearch->ops->snesfunc)(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 = (*linesearch->ops->snesfunc)(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) lambda_update = lambda + PetscRealPart(fty / s);
95 
96     if (PetscIsInfOrNanScalar(lambda_update)) break;
97     if (lambda_update > maxstep) break;
98 
99     /* compute the new state of the line search */
100     lambda_old = lambda;
101     lambda     = lambda_update;
102     fty_old    = fty;
103   }
104   /* construct the solution */
105   ierr = VecCopy(X, W);CHKERRQ(ierr);
106   ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr);
107   if (linesearch->ops->viproject) {
108     ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr);
109   }
110   /* postcheck */
111   ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr);
112   if (changed_y) {
113     ierr = VecAXPY(X, -lambda, Y);CHKERRQ(ierr);
114     if (linesearch->ops->viproject) {
115       ierr = (*linesearch->ops->viproject)(snes, X);CHKERRQ(ierr);
116     }
117   } else {
118     ierr = VecCopy(W, X);CHKERRQ(ierr);
119   }
120   ierr = (*linesearch->ops->snesfunc)(snes,X,F);CHKERRQ(ierr);
121   ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr);
122   if (domainerror) {
123     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
124     PetscFunctionReturn(0);
125   }
126 
127   ierr = SNESLineSearchComputeNorms(linesearch);CHKERRQ(ierr);
128   ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr);
129 
130   ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr);
131 
132   if (monitor) {
133     ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
134     ierr = PetscViewerASCIIPrintf(monitor,"    Line search terminated: lambda = %g, fnorms = %g\n", (double)lambda, (double)gnorm);CHKERRQ(ierr);
135     ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
136   }
137   if (lambda <= steptol) {
138     ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr);
139   }
140   PetscFunctionReturn(0);
141 }
142 
143 #undef __FUNCT__
144 #define __FUNCT__ "SNESLineSearchCreate_CP"
145 /*MC
146    SNESLINESEARCHCP - Critical point line search. This line search assumes that there exists some
147    artificial G(x) for which the SNESFunction F(x) = grad G(x).  Therefore, this line search seeks
148    to find roots of dot(F, Y) via a secant method.
149 
150    Options Database Keys:
151 +  -snes_linesearch_minlambda - the minimum acceptable lambda
152 .  -snes_linesearch_damping - initial trial step length
153 -  -snes_linesearch_max_it  - the maximum number of secant steps performed.
154 
155    Notes:
156    This method is the preferred line search for SNESQN and SNESNCG.
157 
158    Level: advanced
159 
160 .keywords: SNES, SNESLineSearch, damping
161 
162 .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
163 M*/
164 PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch)
165 {
166   PetscFunctionBegin;
167   linesearch->ops->apply          = SNESLineSearchApply_CP;
168   linesearch->ops->destroy        = NULL;
169   linesearch->ops->setfromoptions = NULL;
170   linesearch->ops->reset          = NULL;
171   linesearch->ops->view           = NULL;
172   linesearch->ops->setup          = NULL;
173   linesearch->order               = SNES_LINESEARCH_ORDER_LINEAR;
174 
175   linesearch->max_its = 1;
176   PetscFunctionReturn(0);
177 }
178