Actual source code: linesearchcp.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  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, &ltol, &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) {
 31:     if (monitor) {
 32:       PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
 33:       PetscViewerASCIIPrintf(monitor,"    Line search terminated ended at initial point because dot(F,Y) = %g < atol = %g\n",(double)PetscAbsScalar(fty_old), (double)atol);
 34:       PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
 35:     }
 36:     return(0);
 37:   }

 39:   fty_init = fty_old;

 41:   for (i = 0; i < max_its; i++) {
 42:     /* compute the norm at lambda */
 43:     VecCopy(X, W);
 44:     VecAXPY(W, -lambda, Y);
 45:     if (linesearch->ops->viproject) {
 46:       (*linesearch->ops->viproject)(snes, W);
 47:     }
 48:     (*linesearch->ops->snesfunc)(snes,W,F);
 49:     VecDot(F,Y,&fty);

 51:     delLambda = lambda - lambda_old;

 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 && i > 0) break;
 57:     if (monitor) {
 58:       PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
 59:       PetscViewerASCIIPrintf(monitor,"    Line search: lambdas = [%g, %g], ftys = [%g, %g]\n",(double)lambda, (double)lambda_old, (double)PetscRealPart(fty), (double)PetscRealPart(fty_old));
 60:       PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
 61:     }

 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:       VecCopy(X, W);
 68:       VecAXPY(W, -0.5*(lambda + lambda_old), Y);
 69:       if (linesearch->ops->viproject) {
 70:         (*linesearch->ops->viproject)(snes, W);
 71:       }
 72:       (*linesearch->ops->snesfunc)(snes,W,F);
 73:       VecDot(F, Y, &fty_mid1);
 74:       s    = (3.*fty - 4.*fty_mid1 + fty_old) / delLambda;
 75:     } else {
 76:       VecCopy(X, W);
 77:       VecAXPY(W, -0.5*(lambda + lambda_old), Y);
 78:       if (linesearch->ops->viproject) {
 79:         (*linesearch->ops->viproject)(snes, W);
 80:       }
 81:       (*linesearch->ops->snesfunc)(snes,W,F);
 82:       VecDot(F, Y, &fty_mid1);
 83:       VecCopy(X, W);
 84:       VecAXPY(W, -(lambda + 0.5*(lambda - lambda_old)), Y);
 85:       if (linesearch->ops->viproject) {
 86:         (*linesearch->ops->viproject)(snes, W);
 87:       }
 88:       (*linesearch->ops->snesfunc)(snes, W, F);
 89:       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);

 97:     /* switch directions if we stepped out of bounds */
 98:     if (lambda_update < steptol) lambda_update = lambda + PetscRealPart(fty / s);

100:     if (PetscIsInfOrNanReal(lambda_update)) break;
101:     if (lambda_update > maxstep) break;

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:   VecCopy(X, W);
110:   VecAXPY(W, -lambda, Y);
111:   if (linesearch->ops->viproject) {
112:     (*linesearch->ops->viproject)(snes, W);
113:   }
114:   /* postcheck */
115:   SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);
116:   if (changed_y) {
117:     VecAXPY(X, -lambda, Y);
118:     if (linesearch->ops->viproject) {
119:       (*linesearch->ops->viproject)(snes, X);
120:     }
121:   } else {
122:     VecCopy(W, X);
123:   }
124:   (*linesearch->ops->snesfunc)(snes,X,F);

126:   SNESLineSearchComputeNorms(linesearch);
127:   SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);

129:   SNESLineSearchSetLambda(linesearch, lambda);

131:   if (monitor) {
132:     PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
133:     PetscViewerASCIIPrintf(monitor,"    Line search terminated: lambda = %g, fnorms = %g\n", (double)lambda, (double)gnorm);
134:     PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
135:   }
136:   if (lambda <= steptol) {
137:     SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);
138:   }
139:   return(0);
140: }

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.

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.

153:    Notes:
154:    This method does NOT use the objective function if it is provided with SNESSetObjective().

156:    This method is the preferred line search for SNESQN and SNESNCG.

158:    Level: advanced

160: .keywords: SNES, SNESLineSearch, damping

162: .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
163: M*/
164: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch)
165: {
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;

175:   linesearch->max_its = 1;
176:   return(0);
177: }