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