Actual source code: linesearchbt.c

petsc-3.3-p1 2012-06-15
  1: #include <petsc-private/linesearchimpl.h> /*I  "petscsnes.h"  I*/
  2: #include <petsc-private/snesimpl.h>

  4: typedef struct {
  5:   PetscReal        alpha; /* sufficient decrease parameter */
  6: } SNESLineSearch_BT;

 10: /*@
 11:    SNESLineSearchBTSetAlpha - Sets the descent parameter, alpha, in the BT linesearch variant.

 13:    Input Parameters:
 14: +  linesearch - linesearch context
 15: -  alpha - The descent parameter

 17:    Level: intermediate

 19: .seealso: SNESLineSearchSetLambda(), SNESLineSearchGetTolerances() SNESLINESEARCHBT
 20: @*/
 21: PetscErrorCode SNESLineSearchBTSetAlpha(SNESLineSearch linesearch, PetscReal alpha)
 22: {
 23:   SNESLineSearch_BT  *bt;
 26:   bt = (SNESLineSearch_BT *)linesearch->data;
 27:   bt->alpha = alpha;
 28:   return(0);
 29: }


 34: /*@
 35:    SNESLineSearchBTGetAlpha - Gets the descent parameter, alpha, in the BT linesearch variant.

 37:    Input Parameters:
 38: .  linesearch - linesearch context

 40:    Output Parameters:
 41: .  alpha - The descent parameter

 43:    Level: intermediate

 45: .seealso: SNESLineSearchGetLambda(), SNESLineSearchGetTolerances() SNESLINESEARCHBT
 46: @*/
 47: PetscErrorCode SNESLineSearchBTGetAlpha(SNESLineSearch linesearch, PetscReal *alpha)
 48: {
 49:   SNESLineSearch_BT  *bt;
 52:   bt = (SNESLineSearch_BT *)linesearch->data;
 53:   *alpha = bt->alpha;
 54:   return(0);
 55: }

 59: static PetscErrorCode  SNESLineSearchApply_BT(SNESLineSearch linesearch)
 60: {
 61:   PetscBool      changed_y,changed_w;
 63:   Vec            X,F,Y,W,G;
 64:   SNES           snes;
 65:   PetscReal      fnorm, xnorm, ynorm, gnorm, gnormprev;
 66:   PetscReal      lambda,lambdatemp,lambdaprev,minlambda,maxstep,initslope,alpha,stol;
 67:   PetscReal      t1,t2,a,b,d;
 68: #if defined(PETSC_USE_COMPLEX)
 69:   PetscScalar    cinitslope;
 70: #endif
 71:   PetscBool      domainerror;
 72:   PetscViewer    monitor;
 73:   PetscInt       max_its,count;
 74:   SNESLineSearch_BT  *bt;
 75:   Mat            jac;



 80:   SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);
 81:   SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);
 82:   SNESLineSearchGetLambda(linesearch, &lambda);
 83:   SNESLineSearchGetSNES(linesearch, &snes);
 84:   SNESLineSearchGetMonitor(linesearch, &monitor);
 85:   SNESLineSearchGetTolerances(linesearch,&minlambda,&maxstep,PETSC_NULL,PETSC_NULL,PETSC_NULL,&max_its);
 86:   SNESGetTolerances(snes,PETSC_NULL,PETSC_NULL,&stol,PETSC_NULL,PETSC_NULL);
 87:   bt = (SNESLineSearch_BT *)linesearch->data;

 89:   alpha = bt->alpha;

 91:   SNESGetJacobian(snes, &jac, PETSC_NULL, PETSC_NULL, PETSC_NULL);
 92:   if (!jac) {
 93:     SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix");
 94:   }
 95:   /* precheck */
 96:   SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);
 97:   SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);

 99:   VecNormBegin(Y, NORM_2, &ynorm);
100:   VecNormBegin(X, NORM_2, &xnorm);
101:   VecNormEnd(Y, NORM_2, &ynorm);
102:   VecNormEnd(X, NORM_2, &xnorm);

104:   if (ynorm == 0.0) {
105:     if (monitor) {
106:       PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
107:       PetscViewerASCIIPrintf(monitor,"    Line search: Initial direction and size is 0\n");
108:       PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
109:     }
110:     VecCopy(X,W);
111:     VecCopy(F,G);
112:     SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);
113:     return(0);
114:   }
115:   if (ynorm > maxstep) {        /* Step too big, so scale back */
116:     if (monitor) {
117:       PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
118:       PetscViewerASCIIPrintf(monitor,"    Line search: Scaling step by %14.12e old ynorm %14.12e\n", (maxstep/ynorm),ynorm);
119:       PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
120:     }
121:     VecScale(Y,maxstep/(ynorm));
122:     ynorm = maxstep;
123:   }
124:   MatMult(jac,Y,W);
125: #if defined(PETSC_USE_COMPLEX)
126:   VecDot(F,W,&cinitslope);
127:   initslope = PetscRealPart(cinitslope);
128: #else
129:   VecDot(F,W,&initslope);
130: #endif
131:   if (initslope > 0.0)  initslope = -initslope;
132:   if (initslope == 0.0) initslope = -1.0;

134:   VecWAXPY(W,-lambda,Y,X);
135:   if (linesearch->ops->viproject) {
136:     (*linesearch->ops->viproject)(snes, W);
137:   }
138:   if (snes->nfuncs >= snes->max_funcs) {
139:     PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");
140:     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
141:     SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);
142:     return(0);
143:   }
144:   SNESComputeFunction(snes,W,G);
145:   SNESGetFunctionDomainError(snes, &domainerror);
146:   if (domainerror) {
147:     SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);
148:     return(0);
149:   }
150:   if (linesearch->ops->vinorm) {
151:     gnorm = fnorm;
152:     (*linesearch->ops->vinorm)(snes, G, W, &gnorm);
153:   } else {
154:     VecNorm(G,NORM_2,&gnorm);
155:   }

157:   if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
158:   PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);
159:   if (.5*gnorm*gnorm <= .5*fnorm*fnorm + lambda*alpha*initslope || ynorm < stol*xnorm) { /* Sufficient reduction or step tolerance convergence */
160:     if (monitor) {
161:       PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
162:       PetscViewerASCIIPrintf(monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", fnorm, gnorm);
163:       PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
164:     }
165:   } else {
166:     /* Fit points with quadratic */
167:     lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*lambda*initslope);
168:     lambdaprev = lambda;
169:     gnormprev  = gnorm;
170:     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
171:     if (lambdatemp <= .1*lambda) lambda = .1*lambda;
172:     else                         lambda = lambdatemp;

174:     VecWAXPY(W,-lambda,Y,X);
175:     if (linesearch->ops->viproject) {
176:       (*linesearch->ops->viproject)(snes, W);
177:     }
178:     if (snes->nfuncs >= snes->max_funcs) {
179:       PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);
180:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
181:       SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);
182:       return(0);
183:     }
184:     SNESComputeFunction(snes,W,G);
185:     SNESGetFunctionDomainError(snes, &domainerror);
186:     if (domainerror) {
187:       return(0);
188:     }
189:     if (linesearch->ops->vinorm) {
190:       gnorm = fnorm;
191:       (*linesearch->ops->vinorm)(snes, G, W, &gnorm);
192:     } else {
193:       VecNorm(G,NORM_2,&gnorm);
194:     }
195:     if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
196:     if (monitor) {
197:       PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
198:       PetscViewerASCIIPrintf(monitor,"    Line search: gnorm after quadratic fit %14.12e\n",gnorm);
199:       PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
200:     }
201:     if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
202:       if (monitor) {
203:         PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
204:         PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);
205:         PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
206:       }
207:     } else {
208:       /* Fit points with cubic */
209:       for (count = 0; count < max_its; count++) {
210:         if (lambda <= minlambda) {
211:           if (monitor) {
212:             PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
213:             PetscViewerASCIIPrintf(monitor,"    Line search: unable to find good step length! After %D tries \n",count);
214:             PetscViewerASCIIPrintf(monitor,
215:                                           "    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
216:                                           fnorm, gnorm, ynorm, minlambda, lambda, initslope);
217:             PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
218:           }
219:           SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);
220:           return(0);
221:         }
222:         if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
223:           t1 = .5*(gnorm*gnorm - fnorm*fnorm) - lambda*initslope;
224:           t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
225:           a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
226:           b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
227:           d  = b*b - 3*a*initslope;
228:           if (d < 0.0) d = 0.0;
229:           if (a == 0.0) {
230:             lambdatemp = -initslope/(2.0*b);
231:           } else {
232:             lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a);
233:           }
234:         } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
235:           lambdatemp = -initslope/(gnorm*gnorm - fnorm*fnorm - 2.0*initslope);
236:         } else {
237:           SETERRQ(((PetscObject)linesearch)->comm, PETSC_ERR_SUP, "unsupported line search order for type bt");
238:         }
239:           lambdaprev = lambda;
240:           gnormprev  = gnorm;
241:         if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
242:         if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
243:         else                         lambda     = lambdatemp;
244:         VecWAXPY(W,-lambda,Y,X);
245:         if (snes->nfuncs >= snes->max_funcs) {
246:           PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);
247:           PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",
248:                             fnorm,gnorm,ynorm,lambda,initslope);
249:           SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);
250:           snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
251:           return(0);
252:         }
253:         SNESComputeFunction(snes,W,G);
254:         SNESGetFunctionDomainError(snes, &domainerror);
255:         if (domainerror) {
256:           SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);
257:           return(0);
258:         }
259:         if (linesearch->ops->vinorm) {
260:           gnorm = fnorm;
261:           (*linesearch->ops->vinorm)(snes, G, W, &gnorm);
262:         } else {
263:           VecNorm(G,NORM_2,&gnorm);
264:         }
265:         if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
266:         if (.5*gnorm*gnorm < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
267:           if (monitor) {
268:             PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
269:             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
270:               PetscViewerASCIIPrintf(monitor,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);
271:             } else {
272:               PetscViewerASCIIPrintf(monitor,"    Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",gnorm,lambda);
273:             }
274:             PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
275:           }
276:           break;
277:         } else {
278:           if (monitor) {
279:             PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);
280:             if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
281:               PetscViewerASCIIPrintf(monitor,"    Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);
282:             } else {
283:               PetscViewerASCIIPrintf(monitor,"    Line search: Quadratic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",gnorm,lambda);
284:             }
285:             PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);
286:           }
287:         }
288:       }
289:     }
290:   }

292:   /* postcheck */
293:   SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);
294:   if (changed_y) {
295:     VecWAXPY(W,-lambda,Y,X);
296:     if (linesearch->ops->viproject) {
297:       (*linesearch->ops->viproject)(snes, W);
298:     }
299:   }
300:   if (changed_y || changed_w) { /* recompute the function if the step has changed */
301:     SNESComputeFunction(snes,W,G);
302:     SNESGetFunctionDomainError(snes, &domainerror);
303:     if (domainerror) {
304:       SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);
305:       return(0);
306:     }
307:     if (linesearch->ops->vinorm) {
308:       gnorm = fnorm;
309:       (*linesearch->ops->vinorm)(snes, G, W, &gnorm);
310:     } else {
311:       VecNorm(G,NORM_2,&gnorm);
312:     }
313:     VecNorm(Y,NORM_2,&ynorm);
314:     if (PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");

316:   }

318:   /* copy the solution over */
319:   VecCopy(W, X);
320:   VecCopy(G, F);
321:   VecNorm(X, NORM_2, &xnorm);
322:   SNESLineSearchSetLambda(linesearch, lambda);
323:   SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);
324:   return(0);
325: }

329: PetscErrorCode SNESLineSearchView_BT(SNESLineSearch linesearch, PetscViewer viewer)
330: {
331:   PetscErrorCode    ierr;
332:   PetscBool         iascii;
333:   SNESLineSearch_BT *bt;
335:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
336:   bt = (SNESLineSearch_BT*)linesearch->data;
337:   if (iascii) {
338:     if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) {
339:     PetscViewerASCIIPrintf(viewer, "  interpolation: cubic\n");
340:     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
341:     PetscViewerASCIIPrintf(viewer, "  interpolation: quadratic\n");
342:     }
343:     PetscViewerASCIIPrintf(viewer, "  alpha=%e\n", bt->alpha);
344:   }
345:   return(0);
346: }


351: static PetscErrorCode SNESLineSearchDestroy_BT(SNESLineSearch linesearch)
352: {

356:   PetscFree(linesearch->data);
357:   return(0);
358: }


363: static PetscErrorCode SNESLineSearchSetFromOptions_BT(SNESLineSearch linesearch)
364: {

366:   PetscErrorCode       ierr;
367:   SNESLineSearch_BT    *bt;

370:   bt = (SNESLineSearch_BT*)linesearch->data;

372:   PetscOptionsHead("SNESLineSearch BT options");
373:   PetscOptionsReal("-snes_linesearch_alpha",   "Descent tolerance",        "SNESLineSearchBT", bt->alpha, &bt->alpha, PETSC_NULL);

375:   PetscOptionsTail();
376:   return(0);
377: }


382: /*MC
383:    SNESLINESEARCHBT - Backtracking line search.

385:    This line search finds the minimum of a polynomial fitting of the L2 norm of the
386:    function. If this fit does not satisfy the conditions for progress, the interval shrinks
387:    and the fit is reattempted at most max_it times or until lambda is below minlambda.

389:    Options Database Keys:
390: +  -snes_linesearch_alpha<1e-4> - slope descent parameter
391: .  -snes_linesearch_damping<1.0> - initial step length
392: .  -snes_linesearch_max_it<40> - maximum number of shrinking step
393: .  -snes_linesearch_minlambda<1e-12> - minimum step length allowed
394: -  -snes_linesearch_order<cubic,quadratic> - order of the approximation

396:    Level: advanced

398:    Notes:
399:    This line search is taken from "Numerical Methods for Unconstrained
400:    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.

402: .keywords: SNES, SNESLineSearch, damping

404: .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
405: M*/
407: {

409:   SNESLineSearch_BT  *bt;

413:   linesearch->ops->apply          = SNESLineSearchApply_BT;
414:   linesearch->ops->destroy        = SNESLineSearchDestroy_BT;
415:   linesearch->ops->setfromoptions = SNESLineSearchSetFromOptions_BT;
416:   linesearch->ops->reset          = PETSC_NULL;
417:   linesearch->ops->view           = SNESLineSearchView_BT;
418:   linesearch->ops->setup          = PETSC_NULL;

420:   PetscNewLog(linesearch, SNESLineSearch_BT, &bt);
421:   linesearch->data = (void *)bt;
422:   linesearch->max_its = 40;
423:   linesearch->order = SNES_LINESEARCH_ORDER_CUBIC;
424:   bt->alpha = 1e-4;

426:   return(0);
427: }