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