Actual source code: boundproj.c

  2: #include "boundproj.h"    /*I "tao_solver.h" I*/

  6: static int TaoDestroy_LineSearch(TAO_SOLVER tao,void *linectx)
  7: {
  8:   int  info;
  9:   TAO_LINESEARCH2 *ctx = (TAO_LINESEARCH2 *)linectx;

 11:   TaoFunctionBegin;
 12:   if (ctx->setupcalled==1){
 13:     info=TaoVecDestroy(ctx->W2);CHKERRQ(info);
 14:   }
 15:   info = TaoFree(ctx);CHKERRQ(info);
 16:   TaoFunctionReturn(0);
 17: }

 21: static int TaoSetOptions_LineSearch(TAO_SOLVER tao,void *linectx)
 22: {
 23:   TAO_LINESEARCH2 *ctx = (TAO_LINESEARCH2 *)linectx;
 24:   int            info;

 26:   TaoFunctionBegin;
 27:   info = TaoOptionsHead("Projected line search options");CHKERRQ(info);
 28:   info = TaoOptionInt("-tao_ls_maxfev","max function evals in line search","",ctx->maxfev,&ctx->maxfev,0);CHKERRQ(info);
 29:   info = TaoOptionDouble("-tao_ls_ftol","tol for sufficient decrease","",ctx->ftol,&ctx->ftol,0);CHKERRQ(info);
 30:   info = TaoOptionDouble("-tao_ls_stepmin","lower bound for step","",ctx->stepmin,&ctx->stepmin,0);CHKERRQ(info);
 31:   info = TaoOptionDouble("-tao_ls_stepmax","upper bound for step","",ctx->stepmax,&ctx->stepmax,0);CHKERRQ(info);
 32:   info = TaoOptionsTail();CHKERRQ(info);

 34:   TaoFunctionReturn(0);
 35: }

 39: static int TaoView_LineSearch(TAO_SOLVER tao,void *ctx)
 40: {
 41:   TAO_LINESEARCH2 *ls = (TAO_LINESEARCH2 *)ctx;
 42:   int            info;

 44:   TaoFunctionBegin;
 45:   info = TaoPrintInt(tao,"  Line search: maxf=%d,",ls->maxfev);CHKERRQ(info);
 46:   info = TaoPrintDouble(tao," ftol=%g\n",ls->ftol);CHKERRQ(info);
 47:   TaoFunctionReturn(0);
 48: }

 52: /* @ TaoApply_LineSearch - This routine performs a line search. It
 53:    backtracks until the Armijo conditions are satisfied.

 55:    Input Parameters:
 56: +  tao - TAO_SOLVER context
 57: .  X - current iterate (on output X contains new iterate, X + step*S)
 58: .  S - search direction
 59: .  f - objective function evaluated at X
 60: .  G - gradient evaluated at X
 61: .  W - work vector
 62: .  gdx - inner product of gradient and the direction of the first linear manifold being searched
 63: -  step - initial estimate of step length

 65:    Output parameters:
 66: +  f - objective function evaluated at new iterate, X + step*S
 67: .  G - gradient evaluated at new iterate, X + step*S
 68: .  X - new iterate
 69: -  step - final step length

 71:    Info is set to one of:
 72: .   0 - the line search succeeds; the sufficient decrease
 73:    condition and the directional derivative condition hold

 75:    negative number if an input parameter is invalid
 76: .   -1 -  step < 0 
 77: .   -2 -  ftol < 0 
 78: -   -7 -  maxfev < 0

 80:    positive number > 1 if the line search otherwise terminates
 81: +    2 -  Relative width of the interval of uncertainty is 
 82:          at most rtol.
 83: .    3 -  Maximum number of function evaluations (maxfev) has 
 84:          been reached.
 85: .    4 -  Step is at the lower bound, stepmin.
 86: .    5 -  Step is at the upper bound, stepmax.
 87: .    6 -  Rounding errors may prevent further progress. 
 88:          There may not be a step that satisfies the 
 89:          sufficient decrease and curvature conditions.  
 90:          Tolerances may be too small.
 91: +    7 -  Search direction is not a descent direction.

 93:    Notes:
 94:    This routine is used within the TAO_NLS method.
 95: @ */
 96: static int TaoApply_LineSearch(TAO_SOLVER tao,TaoVec* X,TaoVec* G,TaoVec* S,TaoVec* Gold,double *f,double *f_full,
 97:                         double *step,int *info2,void*ctx)
 98: {
 99:   TAO_LINESEARCH2 *neP = (TAO_LINESEARCH2 *) ctx;
100:   int       info, i;
101:   double zero=0.0;
102:   double finit,gdx,dginit,actred,prered,rho,d1=0,d2=0,d3=0;
103:   TaoVec* XL, *XU, *Xold=neP->W2;
104:   TaoTruth flag;

106:   TaoFunctionBegin;
107:   /* neP->stepmin - lower bound for step */
108:   /* neP->stepmax - upper bound for step */
109:   /* neP->rtol           - relative tolerance for an acceptable step */
110:   /* neP->ftol           - tolerance for sufficient decrease condition */
111:   /* neP->gtol           - tolerance for curvature condition */
112:   /* neP->nfev           - number of function evaluations */
113:   /* neP->maxfev  - maximum number of function evaluations */

115:   /* Check input parameters for errors */

117:   *info2 = 0;
118:     /* After 2 failures, Check that search direction is a descent direction */
119:     /* This test is not really sufficient.  */
120:   if (neP->setupcalled){
121:     info=X->Compatible(neP->W2,&flag); CHKERRQ(info);
122:     if (flag==TAO_FALSE){
123:       info=TaoVecDestroy(neP->W2); CHKERRQ(info);neP->W2=0;
124:       neP->setupcalled=0;
125:     }
126:   }
127:   if (neP->setupcalled==0){
128:     info=X->Clone(&neP->W2); CHKERRQ(info);
129:     Xold=neP->W2;
130:     neP->setupcalled=1;
131:   }
132:   info = TaoGetVariableBounds(tao,&XL,&XU); CHKERRQ(info);
133:   info = Gold->ScaleCopyFrom(-1.0,S); CHKERRQ(info);
134:   info = Gold->BoundGradientProjection(S,XL,X,XU); CHKERRQ(info);
135:   info = Gold->Dot(G,&dginit);CHKERRQ(info);  /* dginit = G^T S */
136:   dginit*=-1;
137:   gdx=dginit;
138:   if (dginit >= zero) {
139:     info = G->CopyFrom(Gold); CHKERRQ(info);
140:     info = X->CopyFrom(Xold); CHKERRQ(info);
141:     info = PetscInfo(tao,"TaoApply_LineSearch2:Search direction not a descent direction\n"); CHKERRQ(info);
142:     *info2 = 7; TaoFunctionReturn(0);
143:   }

145:   info = Xold->CopyFrom(X); CHKERRQ(info);
146:   info = Gold->CopyFrom(G); CHKERRQ(info);
147:   info=TaoGetVariableBounds(tao,&XL,&XU); CHKERRQ(info);
148:   if (*step < zero) {
149:     info = PetscInfo1(tao,"TaoApply_LineSearch:Line search error: step (%g) < 0\n",*step); CHKERRQ(info);
150:     *info2 = -1; TaoFunctionReturn(0);
151:   } else if (neP->ftol < zero) {
152:     info = PetscInfo1(tao,"TaoApply_LineSearch:Line search error: ftol (%g) < 0\n",neP->ftol); CHKERRQ(info);
153:     *info2 = -2; TaoFunctionReturn(0);
154:   } else if (neP->maxfev < zero) {
155:     info = PetscInfo1(tao,"TaoApply_LineSearch:Line search error: maxfev (%d) < 0\n",neP->maxfev); CHKERRQ(info);
156:     *info2 = -7; TaoFunctionReturn(0);
157:   }

159:   /* Initialization */
160:   neP->nfev = 0;
161:   finit = *f;
162:   for (i=0; i< neP->maxfev; i++) {
163:     
164:     /* Force the step to be within the bounds */
165:     *step = TaoMax(*step,neP->stepmin);
166:     *step = TaoMin(*step,neP->stepmax);
167:     
168:     if (0==i){
169:       info = X->Axpy(*step,S);CHKERRQ(info);
170:       info = X->Median(XL,X,XU);CHKERRQ(info);
171:       info = TaoComputeMeritFunctionGradient(tao,X,f,G);
172:       *f_full = *f;

174:       neP->nfev++;
175:       if (info==0 && *f==*f){ /* Successful function evaluation */
176:         actred = *f - finit;
177:         rho = actred/( (*step) * (-TaoAbsScalar(gdx)) );
178:         if (actred < 0 && rho > neP->ftol){
179:           break;
180:         } else{
181:           info=Xold->StepBoundInfo(XL,XU,S,&d3,&d2,&d1);CHKERRQ(info);
182:           info = PetscInfo1(tao,"stepmax = %10f\n",d1); CHKERRQ(info);

184:           *step = TaoMin(*step,d1);
185:           info = Gold->Dot(X,&d1); CHKERRQ(info);
186:           info = Gold->Dot(Xold,&prered); CHKERRQ(info);
187:           prered=d1-prered;
188:           rho = actred/(-TaoAbsScalar(prered));
189:           if (actred < 0 && rho > neP->ftol){
190:             break;
191:           } else{
192:                 *step = (*step)/2;
193:           }
194:         }
195:       } else { /* Function could not be evaluated at this point */
196:         *step = (*step)*0.7;
197:       }

199:     } else {
200:       info = X->Waxpby(*step,S,1.0,Xold);CHKERRQ(info);
201:       info = X->Median(XL,X,XU);CHKERRQ(info);
202:       
203:       info = G->Waxpby(-1,Xold,1.0,X);CHKERRQ(info);
204:       info = G->Dot(Gold,&prered); CHKERRQ(info);
205:       info = TaoComputeMeritFunctionGradient(tao,X,f,G); CHKERRQ(info);
206:       neP->nfev++;
207:       if (info==0 && *f==*f){ /* Successful function evaluation */
208:         actred = *f - finit;
209:         rho = actred/(-TaoAbsScalar(prered));
210:         /* 
211:            If sufficient progress has been obtained, accept the
212:            point.   Prered could be positive or negative.  
213:            Otherwise, backtrack. 
214:         */

216:         if (actred < 0 && rho > neP->ftol){
217:           info = PetscInfo4(tao,"TaoApply_LineSearch: steplength: %g, actual reduction: %g (hopefully positive), Predicted reduction: %g, rho: %g\n",*step,-actred,prered,rho); CHKERRQ(info);
218:           break;
219:         } else {
220:           info = PetscInfo4(tao,"TaoApply_LineSearch: steplength: %g, actual reduction: %g (hopefully positive), Predicted reduction: %g, rho: %g\n",*step,-actred,prered,rho); CHKERRQ(info);
221:           if (i<5){
222:             *step = (*step)/2;
223:           } else {
224:             *step = (*step)/2;
225:           }
226:         }
227:       } else { /* Function could not be evaluated at this point */
228:           info = PetscInfo(tao,"TaoApply_LineSearch: Problem in function evaluation\n"); CHKERRQ(info);
229:         *step = (*step)*0.7;
230:       }
231:     }

233:   }
234:   /* Convergence testing */
235:   
236:   if (*step <= neP->stepmin || *step >= neP->stepmax) {
237:     *info2 = 6;
238:     info = PetscInfo(tao,"TaoApply_LineSearch:Rounding errors may prevent further progress.  May not be a step satisfying\n"); CHKERRQ(info);
239:     info = PetscInfo(tao,"TaoApply_LineSearch:sufficient decrease and curvature conditions. Tolerances may be too small.\n"); CHKERRQ(info);
240:   }
241:   if (*step == neP->stepmax) {
242:     info = PetscInfo1(tao,"TaoApply_LineSearch:Step is at the upper bound, stepmax (%g)\n",neP->stepmax); CHKERRQ(info);
243:     *info2 = 5;
244:   }
245:   if (*step == neP->stepmin) {
246:     info = PetscInfo1(tao,"TaoApply_LineSearch:Step is at the lower bound, stepmin (%g)\n",neP->stepmin); CHKERRQ(info);
247:     *info2 = 4;
248:   }
249:   if (neP->nfev >= neP->maxfev) {
250:     info = PetscInfo2(tao,"TaoApply_LineSearch:Number of line search function evals (%d) > maximum (%d)\n",neP->nfev,neP->maxfev); CHKERRQ(info);
251:     *info2 = 3;
252:   }
253:   if ((neP->bracket) && (neP->stepmax - neP->stepmin <= neP->rtol*neP->stepmax)){
254:     info = PetscInfo1(tao,"TaoApply_LineSearch:Relative width of interval of uncertainty is at most rtol (%g)\n",neP->rtol); CHKERRQ(info);
255:     *info2 = 2;
256:   }
257:   /*
258:   if ((*f <= ftest1) && (TaoAbsDouble(dg) <= neP->gtol*(-dginit))) {
259:     info = PetscLogInfo((tao,"TaoApply_LineSearch:Line search success: Sufficient decrease and directional deriv conditions hold\n")); CHKERRQ(info);
260:     *info2 = 1;
261:   }
262:   */
263:   
264:   /* Finish computations */
265:   /*
266:   info = PetscLogInfo((tao,"TaoApply_LineSearch:%d function evals in line search, step = %10.4f\n",neP->nfev,*step)); CHKERRQ(info);
267:   */
268:   TaoFunctionReturn(0);
269: }

273: /*@C
274:    TaoCreateProjectedLineSearch - Create a line search

276:    Input Parameters:
277: .  tao - TAO_SOLVER context

279:    Note:
280:    This routine is used within the following TAO bound constrained 
281:    minimization solvers: TRON (tao_tron) and limited memory variable metric
282:    (tao_blmvm). This line search projects points onto the feasible, bound
283:    constrained region.  It only enforces the Armijo descent condition.

285:    Level: developer

287: .keywords: TAO_SOLVER, linesearch
288: @*/
289: int TaoCreateProjectedLineSearch(TAO_SOLVER tao)
290: {
291:   int info;
292:   TAO_LINESEARCH2 *neP;

294:   TaoFunctionBegin;

296:   info = TaoNew(TAO_LINESEARCH2,&neP); CHKERRQ(info);
297:   info = PetscLogObjectMemory(tao,sizeof(TAO_LINESEARCH2)); CHKERRQ(info);

299:   neP->ftol                  = 0.001;
300:   neP->rtol                  = 0.0;
301:   neP->gtol                  = 0.0;
302:   neP->stepmin                  = 1.0e-30;
303:   neP->stepmax                  = 1.0e+20;
304:   neP->nfev                  = 0; 
305:   neP->bracket                  = 0; 
306:   neP->infoc              = 1;
307:   neP->maxfev                  = 30;
308:   neP->W2                 = TAO_NULL;
309:   neP->setupcalled        = 0;

311:   info = TaoSetLineSearch(tao,0,
312:                           TaoSetOptions_LineSearch,
313:                           TaoApply_LineSearch,
314:                           TaoView_LineSearch,
315:                           TaoDestroy_LineSearch,
316:                           (void *) neP);CHKERRQ(info);

318:   TaoFunctionReturn(0);
319: }