Actual source code: linesearchnleqerr.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  1: #include <petsc/private/linesearchimpl.h> /*I  "petscsnes.h"  I*/
  2: #include <petsc/private/snesimpl.h>

  4: typedef struct {
  5:   PetscReal norm_delta_x_prev; /* norm of previous update */
  6:   PetscReal norm_bar_delta_x_prev; /* norm of previous bar update */
  7:   PetscReal mu_curr; /* current local Lipschitz estimate */
  8:   PetscReal lambda_prev; /* previous step length: for some reason SNESLineSearchGetLambda returns 1 instead of the previous step length */
  9: } SNESLineSearch_NLEQERR;

 11: static PetscBool NLEQERR_cited = PETSC_FALSE;
 12: static const char NLEQERR_citation[] = "@book{deuflhard2011,\n"
 13:                                "  title = {Newton Methods for Nonlinear Problems},\n"
 14:                                "  author = {Peter Deuflhard},\n"
 15:                                "  volume = 35,\n"
 16:                                "  year = 2011,\n"
 17:                                "  isbn = {978-3-642-23898-7},\n"
 18:                                "  doi  = {10.1007/978-3-642-23899-4},\n"
 19:                                "  publisher = {Springer-Verlag},\n"
 20:                                "  address = {Berlin, Heidelberg}\n}\n";

 24: static PetscErrorCode SNESLineSearchReset_NLEQERR(SNESLineSearch linesearch)
 25: {
 26:   SNESLineSearch_NLEQERR *nleqerr = (SNESLineSearch_NLEQERR*)linesearch->data;

 29:   nleqerr->mu_curr               = 0.0;
 30:   nleqerr->norm_delta_x_prev     = -1.0;
 31:   nleqerr->norm_bar_delta_x_prev = -1.0;
 32:   return(0);
 33: }

 37: static PetscErrorCode  SNESLineSearchApply_NLEQERR(SNESLineSearch linesearch)
 38: {
 39:   PetscBool              changed_y,changed_w;
 40:   PetscErrorCode         ierr;
 41:   Vec                    X,F,Y,W,G;
 42:   SNES                   snes;
 43:   PetscReal              fnorm, xnorm, ynorm, gnorm, wnorm;
 44:   PetscReal              lambda, minlambda, stol;
 45:   PetscViewer            monitor;
 46:   PetscInt               max_its, count, snes_iteration;
 47:   PetscReal              theta, mudash, lambdadash;
 48:   SNESLineSearch_NLEQERR *nleqerr = (SNESLineSearch_NLEQERR*)linesearch->data;
 49:   KSPConvergedReason     kspreason;

 52:   PetscCitationsRegister(NLEQERR_citation, &NLEQERR_cited);

 54:   SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);
 55:   SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);
 56:   SNESLineSearchGetLambda(linesearch, &lambda);
 57:   SNESLineSearchGetSNES(linesearch, &snes);
 58:   SNESLineSearchGetDefaultMonitor(linesearch, &monitor);
 59:   SNESLineSearchGetTolerances(linesearch,&minlambda,NULL,NULL,NULL,NULL,&max_its);
 60:   SNESGetTolerances(snes,NULL,NULL,&stol,NULL,NULL);

 62:   /* reset the state of the Lipschitz estimates */
 63:   SNESGetIterationNumber(snes, &snes_iteration);
 64:   if (!snes_iteration) {
 65:     SNESLineSearchReset_NLEQERR(linesearch);
 66:   }

 68:   /* precheck */
 69:   SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);
 70:   SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED);

 72:   VecNormBegin(Y, NORM_2, &ynorm);
 73:   VecNormBegin(X, NORM_2, &xnorm);
 74:   VecNormEnd(Y, NORM_2, &ynorm);
 75:   VecNormEnd(X, NORM_2, &xnorm);

 77:   /* Note: Y is *minus* the Newton step. For whatever reason PETSc doesn't solve with the minus on  the RHS. */

 79:   if (ynorm == 0.0) {
 80:     if (monitor) {
 81:       PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
 82:       PetscViewerASCIIPrintf(monitor,"    Line search: Initial direction and size is 0\n");
 83:       PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
 84:     }
 85:     VecCopy(X,W);
 86:     VecCopy(F,G);
 87:     SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm);
 88:     SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);
 89:     return(0);
 90:   }

 92:   /* At this point, we've solved the Newton system for delta_x, and we assume that
 93:      its norm is greater than the solution tolerance (otherwise we wouldn't be in
 94:      here). So let's go ahead and estimate the Lipschitz constant. 

 96:      W contains bar_delta_x_prev at this point. */

 98:   if (monitor) {
 99:     PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
100:     PetscViewerASCIIPrintf(monitor,"    Line search: norm of Newton step: %14.12e\n", (double) ynorm);
101:     PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
102:   }

104:   /* this needs information from a previous iteration, so can't do it on the first one */
105:   if (nleqerr->norm_delta_x_prev > 0 && nleqerr->norm_bar_delta_x_prev > 0) {
106:     VecWAXPY(G, +1.0, Y, W); /* bar_delta_x - delta_x; +1 because Y is -delta_x */
107:     VecNormBegin(G, NORM_2, &gnorm);
108:     VecNormEnd(G, NORM_2, &gnorm);

110:     nleqerr->mu_curr = nleqerr->lambda_prev * (nleqerr->norm_delta_x_prev * nleqerr->norm_bar_delta_x_prev) / (gnorm * ynorm);
111:     lambda = PetscMin(1.0, nleqerr->mu_curr);

113:     if (monitor) {
114:       PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
115:       PetscViewerASCIIPrintf(monitor,"    Line search: Lipschitz estimate: %14.12e; lambda: %14.12e\n", (double) nleqerr->mu_curr, (double) lambda);
116:       PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
117:     }
118:   } else {
119:     lambda = linesearch->damping;
120:   }

122:   /* The main while loop of the algorithm. 
123:      At the end of this while loop, G should have the accepted new X in it. */

125:   count = 0;
126:   while (PETSC_TRUE) {
127:     if (monitor) {
128:       PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
129:       PetscViewerASCIIPrintf(monitor,"    Line search: entering iteration with lambda: %14.12e\n", lambda);
130:       PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
131:     }

133:     /* Check that we haven't performed too many iterations */
134:     count += 1;
135:     if (count >= max_its) {
136:       if (monitor) {
137:         PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
138:         PetscViewerASCIIPrintf(monitor,"    Line search: maximum iterations reached\n");
139:         PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
140:       }
141:       SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);
142:       return(0);
143:     }

145:     /* Now comes the Regularity Test. */
146:     if (lambda <= minlambda) {
147:       /* This isn't what is suggested by Deuflhard, but it works better in my experience */
148:       if (monitor) {
149:         PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
150:         PetscViewerASCIIPrintf(monitor,"    Line search: lambda has reached lambdamin, taking full Newton step\n");
151:         PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
152:       }
153:       lambda = 1.0;
154:       VecWAXPY(G, -lambda, Y, X);

156:       /* and clean up the state for next time */
157:       SNESLineSearchReset_NLEQERR(linesearch);
158:       /*
159:          The clang static analyzer detected a problem here; once the loop is broken the values
160:          nleqerr->norm_delta_x_prev     = ynorm;
161:          nleqerr->norm_bar_delta_x_prev = wnorm;
162:          are set, but wnorm has not even been computed.
163:          I don't know if this is the correct fix but by setting ynorm and wnorm to -1.0 at
164:          least the linesearch object is kept in the state set by the SNESLineSearchReset_NLEQERR() call above
165:       */
166:       ynorm = wnorm = -1.0;
167:       break;
168:     }

170:     /* Compute new trial iterate */
171:     VecWAXPY(W, -lambda, Y, X);
172:     SNESComputeFunction(snes, W, G);

174:     /* Solve linear system for bar_delta_x_curr: old Jacobian, new RHS. Note absence of minus sign, compared to Deuflhard, in keeping with PETSc convention */
175:     KSPSolve(snes->ksp, G, W);
176:     KSPGetConvergedReason(snes->ksp, &kspreason);
177:     if (kspreason < 0) {
178:       PetscInfo(snes,"Solution for \\bar{delta x}^{k+1} failed.");
179:     }

181:     /* W now contains -bar_delta_x_curr. */

183:     VecNorm(W, NORM_2, &wnorm);
184:     if (monitor) {
185:       PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
186:       PetscViewerASCIIPrintf(monitor,"    Line search: norm of simplified Newton update: %14.12e\n", (double) wnorm);
187:       PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
188:     }

190:     /* compute the monitoring quantities theta and mudash. */

192:     theta = wnorm / ynorm;

194:     VecWAXPY(G, -(1.0 - lambda), Y, W);
195:     VecNorm(G, NORM_2, &gnorm);

197:     mudash = (0.5 * ynorm * lambda * lambda) / gnorm;

199:     /* Check for termination of the linesearch */
200:     if (theta >= 1.0) {
201:       /* need to go around again with smaller lambda */
202:       if (monitor) {
203:         PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
204:         PetscViewerASCIIPrintf(monitor,"    Line search: monotonicity check failed, ratio: %14.12e\n", (double) theta);
205:         PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
206:       }
207:       lambda = PetscMin(mudash, 0.5 * lambda);
208:       lambda = PetscMax(lambda, minlambda);
209:       /* continue through the loop, i.e. go back to regularity test */
210:     } else {
211:       /* linesearch terminated */
212:       lambdadash = PetscMin(1.0, mudash);

214:       if (lambdadash == 1.0 && lambda == 1.0 && wnorm <= stol) {
215:         /* store the updated state, X - Y - W, in G:
216:            I need to keep W for the next linesearch */
217:         VecCopy(X, G);
218:         VecAXPY(G, -1.0, Y);
219:         VecAXPY(G, -1.0, W);
220:         break;
221:       }

223:       /* Deuflhard suggests to add the following:
224:       else if (lambdadash >= 4.0 * lambda) {
225:         lambda = lambdadash;
226:       }
227:       to continue through the loop, i.e. go back to regularity test.
228:       I deliberately exclude this, as I have practical experience of this
229:       getting stuck in infinite loops (on e.g. an Allen--Cahn problem). */

231:       else {
232:         /* accept iterate without adding on, i.e. don't use bar_delta_x;
233:            again, I need to keep W for the next linesearch */
234:         VecWAXPY(G, -lambda, Y, X);
235:         break;
236:       }
237:     }
238:   }

240:   if (linesearch->ops->viproject) {
241:     (*linesearch->ops->viproject)(snes, G);
242:   }

244:   /* W currently contains -bar_delta_u. Scale it so that it contains bar_delta_u. */
245:   VecScale(W, -1.0);

247:   /* postcheck */
248:   SNESLineSearchPostCheck(linesearch,X,Y,G,&changed_y,&changed_w);
249:   if (changed_y || changed_w) {
250:     SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_USER);
251:     PetscInfo(snes,"Changing the search direction here doesn't make sense.\n");
252:     return(0);
253:   }

255:   /* copy the solution and information from this iteration over */
256:   nleqerr->norm_delta_x_prev     = ynorm;
257:   nleqerr->norm_bar_delta_x_prev = wnorm;
258:   nleqerr->lambda_prev           = lambda;

260:   VecCopy(G, X);
261:   SNESComputeFunction(snes, X, F);
262:   VecNorm(X, NORM_2, &xnorm);
263:   VecNorm(F, NORM_2, &fnorm);
264:   SNESLineSearchSetLambda(linesearch, lambda);
265:   SNESLineSearchSetNorms(linesearch, xnorm, fnorm, ynorm);
266:   return(0);
267: }

271: PetscErrorCode SNESLineSearchView_NLEQERR(SNESLineSearch linesearch, PetscViewer viewer)
272: {
273:   PetscErrorCode          ierr;
274:   PetscBool               iascii;
275:   SNESLineSearch_NLEQERR *nleqerr;

278:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
279:   nleqerr   = (SNESLineSearch_NLEQERR*)linesearch->data;
280:   if (iascii) {
281:     PetscViewerASCIIPrintf(viewer, "  NLEQ-ERR affine-covariant linesearch");
282:     PetscViewerASCIIPrintf(viewer, "  current local Lipschitz estimate omega=%e\n", (double)nleqerr->mu_curr);
283:   }
284:   return(0);
285: }

289: static PetscErrorCode SNESLineSearchDestroy_NLEQERR(SNESLineSearch linesearch)
290: {

294:   PetscFree(linesearch->data);
295:   return(0);
296: }

300: /*MC
301:    SNESLINESEARCHNLEQERR - Error-oriented affine-covariant globalised Newton algorithm of Deuflhard (2011).

303:    This linesearch is intended for Newton-type methods which are affine covariant. Affine covariance
304:    means that Newton's method will give the same iterations for F(x) = 0 and AF(x) = 0 for a nonsingular
305:    matrix A. This is a fundamental property; the philosophy of this linesearch is that globalisations
306:    of Newton's method should carefully preserve it.

308:    For a discussion of the theory behind this algorithm, see

310:    @book{deuflhard2011,
311:      title={Newton Methods for Nonlinear Problems},
312:      author={Deuflhard, P.},
313:      volume={35},
314:      year={2011},
315:      publisher={Springer-Verlag},
316:      address={Berlin, Heidelberg}
317:    }

319:    Pseudocode is given on page 148.

321:    Options Database Keys:
322: +  -snes_linesearch_damping<1.0> - initial step length
323: -  -snes_linesearch_minlambda<1e-12> - minimum step length allowed

325:    Contributed by Patrick Farrell <patrick.farrell@maths.ox.ac.uk>

327:    Level: advanced

329: .keywords: SNES, SNESLineSearch, damping

331: .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
332: M*/
333: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NLEQERR(SNESLineSearch linesearch)
334: {

336:   SNESLineSearch_NLEQERR *nleqerr;
337:   PetscErrorCode    ierr;

340:   linesearch->ops->apply          = SNESLineSearchApply_NLEQERR;
341:   linesearch->ops->destroy        = SNESLineSearchDestroy_NLEQERR;
342:   linesearch->ops->setfromoptions = NULL;
343:   linesearch->ops->reset          = SNESLineSearchReset_NLEQERR;
344:   linesearch->ops->view           = SNESLineSearchView_NLEQERR;
345:   linesearch->ops->setup          = NULL;

347:   PetscNewLog(linesearch,&nleqerr);

349:   linesearch->data    = (void*)nleqerr;
350:   linesearch->max_its = 40;
351:   SNESLineSearchReset_NLEQERR(linesearch);
352:   return(0);
353: }