Actual source code: tron.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: #include <../src/tao/bound/impls/tron/tron.h>
  2: #include <petsc-private/kspimpl.h>
  3: #include <petsc-private/matimpl.h>
  4: #include <../src/tao/matrix/submatfree.h>


  7: /* TRON Routines */
  8: static PetscErrorCode TronGradientProjections(Tao,TAO_TRON*);
  9: /*------------------------------------------------------------*/
 12: static PetscErrorCode TaoDestroy_TRON(Tao tao)
 13: {
 14:   TAO_TRON       *tron = (TAO_TRON *)tao->data;

 18:   VecDestroy(&tron->X_New);
 19:   VecDestroy(&tron->G_New);
 20:   VecDestroy(&tron->Work);
 21:   VecDestroy(&tron->DXFree);
 22:   VecDestroy(&tron->R);
 23:   VecDestroy(&tron->diag);
 24:   VecScatterDestroy(&tron->scatter);
 25:   ISDestroy(&tron->Free_Local);
 26:   MatDestroy(&tron->H_sub);
 27:   MatDestroy(&tron->Hpre_sub);
 28:   PetscFree(tao->data);
 29:   return(0);
 30: }

 32: /*------------------------------------------------------------*/
 35: static PetscErrorCode TaoSetFromOptions_TRON(Tao tao)
 36: {
 37:   TAO_TRON       *tron = (TAO_TRON *)tao->data;
 39:   PetscBool      flg;

 42:   PetscOptionsHead("Newton Trust Region Method for bound constrained optimization");
 43:   PetscOptionsInt("-tao_tron_maxgpits","maximum number of gradient projections per TRON iterate","TaoSetMaxGPIts",tron->maxgpits,&tron->maxgpits,&flg);
 44:   PetscOptionsTail();
 45:   TaoLineSearchSetFromOptions(tao->linesearch);
 46:   KSPSetFromOptions(tao->ksp);
 47:   return(0);
 48: }

 50: /*------------------------------------------------------------*/
 53: static PetscErrorCode TaoView_TRON(Tao tao, PetscViewer viewer)
 54: {
 55:   TAO_TRON         *tron = (TAO_TRON *)tao->data;
 56:   PetscBool        isascii;
 57:   PetscErrorCode   ierr;

 60:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 61:   if (isascii) {
 62:     PetscViewerASCIIPushTab(viewer);
 63:     PetscViewerASCIIPrintf(viewer,"Total PG its: %D,",tron->total_gp_its);
 64:     PetscViewerASCIIPrintf(viewer,"PG tolerance: %g \n",(double)tron->pg_ftol);
 65:     PetscViewerASCIIPopTab(viewer);
 66:   }
 67:   return(0);
 68: }


 71: /* ---------------------------------------------------------- */
 74: static PetscErrorCode TaoSetup_TRON(Tao tao)
 75: {
 77:   TAO_TRON       *tron = (TAO_TRON *)tao->data;


 81:   /* Allocate some arrays */
 82:   VecDuplicate(tao->solution, &tron->diag);
 83:   VecDuplicate(tao->solution, &tron->X_New);
 84:   VecDuplicate(tao->solution, &tron->G_New);
 85:   VecDuplicate(tao->solution, &tron->Work);
 86:   VecDuplicate(tao->solution, &tao->gradient);
 87:   VecDuplicate(tao->solution, &tao->stepdirection);
 88:   if (!tao->XL) {
 89:       VecDuplicate(tao->solution, &tao->XL);
 90:       VecSet(tao->XL, PETSC_NINFINITY);
 91:   }
 92:   if (!tao->XU) {
 93:       VecDuplicate(tao->solution, &tao->XU);
 94:       VecSet(tao->XU, PETSC_INFINITY);
 95:   }
 96:   return(0);
 97: }



103: static PetscErrorCode TaoSolve_TRON(Tao tao)
104: {
105:   TAO_TRON                     *tron = (TAO_TRON *)tao->data;
106:   PetscErrorCode               ierr;
107:   PetscInt                     iter=0,its;
108:   TaoConvergedReason           reason = TAO_CONTINUE_ITERATING;
109:   TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING;
110:   PetscReal                    prered,actred,delta,f,f_new,rhok,gdx,xdiff,stepsize;

113:   tron->pgstepsize=1.0;
114:   tao->trust = tao->trust0;
115:   /*   Project the current point onto the feasible set */
116:   TaoComputeVariableBounds(tao);
117:   VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);
118:   TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);

120:   TaoComputeObjectiveAndGradient(tao,tao->solution,&tron->f,tao->gradient);
121:   ISDestroy(&tron->Free_Local);

123:   VecWhichBetween(tao->XL,tao->solution,tao->XU,&tron->Free_Local);

125:   /* Project the gradient and calculate the norm */
126:   VecBoundGradientProjection(tao->gradient,tao->solution, tao->XL, tao->XU, tao->gradient);
127:   VecNorm(tao->gradient,NORM_2,&tron->gnorm);

129:   if (PetscIsInfOrNanReal(tron->f) || PetscIsInfOrNanReal(tron->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf pr NaN");
130:   if (tao->trust <= 0) {
131:     tao->trust=PetscMax(tron->gnorm*tron->gnorm,1.0);
132:   }

134:   tron->stepsize=tao->trust;
135:   TaoMonitor(tao, iter, tron->f, tron->gnorm, 0.0, tron->stepsize, &reason);
136:   while (reason==TAO_CONTINUE_ITERATING){

138:     TronGradientProjections(tao,tron);
139:     f=tron->f; delta=tao->trust;
140:     tron->n_free_last = tron->n_free;
141:     TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);

143:     ISGetSize(tron->Free_Local, &tron->n_free);

145:     /* If no free variables */
146:     if (tron->n_free == 0) {
147:       actred=0;
148:       PetscInfo(tao,"No free variables in tron iteration.");
149:       break;

151:     }
152:     /* use free_local to mask/submat gradient, hessian, stepdirection */
153:     TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->R);
154:     TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->DXFree);
155:     VecSet(tron->DXFree,0.0);
156:     VecScale(tron->R, -1.0);
157:     TaoMatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub);
158:     if (tao->hessian == tao->hessian_pre) {
159:       MatDestroy(&tron->Hpre_sub);
160:       PetscObjectReference((PetscObject)(tron->H_sub));
161:       tron->Hpre_sub = tron->H_sub;
162:     } else {
163:       TaoMatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type,&tron->Hpre_sub);
164:     }
165:     KSPReset(tao->ksp);
166:     KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub);
167:     while (1) {

169:       /* Approximately solve the reduced linear system */
170:       KSPSTCGSetRadius(tao->ksp,delta);

172:       KSPSolve(tao->ksp, tron->R, tron->DXFree);
173:       KSPGetIterationNumber(tao->ksp,&its);
174:       tao->ksp_its+=its;
175:       VecSet(tao->stepdirection,0.0);

177:       /* Add dxfree matrix to compute step direction vector */
178:       VecISAXPY(tao->stepdirection,tron->Free_Local,1.0,tron->DXFree);
179:       if (0) {
180:         PetscReal rhs,stepnorm;
181:         VecNorm(tron->R,NORM_2,&rhs);
182:         VecNorm(tron->DXFree,NORM_2,&stepnorm);
183:         PetscPrintf(PETSC_COMM_WORLD,"|rhs|=%g\t|s|=%g\n",(double)rhs,(double)stepnorm);
184:       }


187:       VecDot(tao->gradient, tao->stepdirection, &gdx);
188:       PetscInfo1(tao,"Expected decrease in function value: %14.12e\n",(double)gdx);

190:       VecCopy(tao->solution, tron->X_New);
191:       VecCopy(tao->gradient, tron->G_New);

193:       stepsize=1.0;f_new=f;

195:       TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
196:       TaoLineSearchApply(tao->linesearch, tron->X_New, &f_new, tron->G_New, tao->stepdirection,&stepsize,&ls_reason);
197:       TaoAddLineSearchCounts(tao);

199:       MatMult(tao->hessian, tao->stepdirection, tron->Work);
200:       VecAYPX(tron->Work, 0.5, tao->gradient);
201:       VecDot(tao->stepdirection, tron->Work, &prered);
202:       actred = f_new - f;
203:       if (actred<0) {
204:         rhok=PetscAbs(-actred/prered);
205:       } else {
206:         rhok=0.0;
207:       }

209:       /* Compare actual improvement to the quadratic model */
210:       if (rhok > tron->eta1) { /* Accept the point */
211:         /* d = x_new - x */
212:         VecCopy(tron->X_New, tao->stepdirection);
213:         VecAXPY(tao->stepdirection, -1.0, tao->solution);

215:         VecNorm(tao->stepdirection, NORM_2, &xdiff);
216:         xdiff *= stepsize;

218:         /* Adjust trust region size */
219:         if (rhok < tron->eta2 ){
220:           delta = PetscMin(xdiff,delta)*tron->sigma1;
221:         } else if (rhok > tron->eta4 ){
222:           delta= PetscMin(xdiff,delta)*tron->sigma3;
223:         } else if (rhok > tron->eta3 ){
224:           delta=PetscMin(xdiff,delta)*tron->sigma2;
225:         }
226:         VecBoundGradientProjection(tron->G_New,tron->X_New, tao->XL, tao->XU, tao->gradient);
227:         if (tron->Free_Local) {
228:           ISDestroy(&tron->Free_Local);
229:         }
230:         VecWhichBetween(tao->XL, tron->X_New, tao->XU, &tron->Free_Local);
231:         f=f_new;
232:         VecNorm(tao->gradient,NORM_2,&tron->gnorm);
233:         VecCopy(tron->X_New, tao->solution);
234:         VecCopy(tron->G_New, tao->gradient);
235:         break;
236:       }
237:       else if (delta <= 1e-30) {
238:         break;
239:       }
240:       else {
241:         delta /= 4.0;
242:       }
243:     } /* end linear solve loop */


246:     tron->f=f; tron->actred=actred; tao->trust=delta;
247:     iter++;
248:     TaoMonitor(tao, iter, tron->f, tron->gnorm, 0.0, delta, &reason);
249:   }  /* END MAIN LOOP  */

251:   return(0);
252: }


257: static PetscErrorCode TronGradientProjections(Tao tao,TAO_TRON *tron)
258: {
259:   PetscErrorCode                 ierr;
260:   PetscInt                       i;
261:   TaoLineSearchConvergedReason ls_reason;
262:   PetscReal                      actred=-1.0,actred_max=0.0;
263:   PetscReal                      f_new;
264:   /*
265:      The gradient and function value passed into and out of this
266:      routine should be current and correct.

268:      The free, active, and binding variables should be already identified
269:   */
271:   if (tron->Free_Local) {
272:     ISDestroy(&tron->Free_Local);
273:   }
274:   VecWhichBetween(tao->XL,tao->solution,tao->XU,&tron->Free_Local);

276:   for (i=0;i<tron->maxgpits;i++){

278:     if ( -actred <= (tron->pg_ftol)*actred_max) break;

280:     tron->gp_iterates++; tron->total_gp_its++;
281:     f_new=tron->f;

283:     VecCopy(tao->gradient, tao->stepdirection);
284:     VecScale(tao->stepdirection, -1.0);
285:     TaoLineSearchSetInitialStepLength(tao->linesearch,tron->pgstepsize);
286:     TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection,
287:                               &tron->pgstepsize, &ls_reason);
288:     TaoAddLineSearchCounts(tao);


291:     /* Update the iterate */
292:     actred = f_new - tron->f;
293:     actred_max = PetscMax(actred_max,-(f_new - tron->f));
294:     tron->f = f_new;
295:     if (tron->Free_Local) {
296:       ISDestroy(&tron->Free_Local);
297:     }
298:     VecWhichBetween(tao->XL,tao->solution,tao->XU,&tron->Free_Local);
299:   }

301:   return(0);
302: }

306: static PetscErrorCode TaoComputeDual_TRON(Tao tao, Vec DXL, Vec DXU) {

308:   TAO_TRON       *tron = (TAO_TRON *)tao->data;

315:   if (!tron->Work || !tao->gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");

317:   VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tron->Work);
318:   VecCopy(tron->Work,DXL);
319:   VecAXPY(DXL,-1.0,tao->gradient);
320:   VecSet(DXU,0.0);
321:   VecPointwiseMax(DXL,DXL,DXU);

323:   VecCopy(tao->gradient,DXU);
324:   VecAXPY(DXU,-1.0,tron->Work);
325:   VecSet(tron->Work,0.0);
326:   VecPointwiseMin(DXU,tron->Work,DXU);
327:   return(0);
328: }

330: /*------------------------------------------------------------*/
331: /*MC
332:   TAOTRON - The TRON algorithm is an active-set Newton trust region method
333:   for bound-constrained minimization.

335:   Options Database Keys:
336: + -tao_tron_maxgpits - maximum number of gradient projections per TRON iterate
337: - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets

339:   Level: beginner
340: M*/
341: EXTERN_C_BEGIN
344: PetscErrorCode TaoCreate_TRON(Tao tao)
345: {
346:   TAO_TRON       *tron;
348:   const char     *morethuente_type = TAOLINESEARCHMT;

351:   tao->ops->setup = TaoSetup_TRON;
352:   tao->ops->solve = TaoSolve_TRON;
353:   tao->ops->view = TaoView_TRON;
354:   tao->ops->setfromoptions = TaoSetFromOptions_TRON;
355:   tao->ops->destroy = TaoDestroy_TRON;
356:   tao->ops->computedual = TaoComputeDual_TRON;

358:   PetscNewLog(tao,&tron);

360:   tao->max_it = 50;
361: #if defined(PETSC_USE_REAL_SINGLE)
362:   tao->fatol = 1e-5;
363:   tao->frtol = 1e-5;
364:   tao->steptol = 1e-6;
365: #else
366:   tao->fatol = 1e-10;
367:   tao->frtol = 1e-10;
368:   tao->steptol = 1e-12;
369: #endif
370:   tao->data = (void*)tron;
371:   tao->trust0       = 1.0;

373:   /* Initialize pointers and variables */
374:   tron->n            = 0;
375:   tron->maxgpits     = 3;
376:   tron->pg_ftol      = 0.001;

378:   tron->eta1         = 1.0e-4;
379:   tron->eta2         = 0.25;
380:   tron->eta3         = 0.50;
381:   tron->eta4         = 0.90;

383:   tron->sigma1       = 0.5;
384:   tron->sigma2       = 2.0;
385:   tron->sigma3       = 4.0;

387:   tron->gp_iterates  = 0; /* Cumulative number */
388:   tron->total_gp_its = 0;
389:   tron->n_free       = 0;

391:   tron->DXFree=NULL;
392:   tron->R=NULL;
393:   tron->X_New=NULL;
394:   tron->G_New=NULL;
395:   tron->Work=NULL;
396:   tron->Free_Local=NULL;
397:   tron->H_sub=NULL;
398:   tron->Hpre_sub=NULL;
399:   tao->subset_type = TAO_SUBSET_SUBVEC;

401:   TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);
402:   TaoLineSearchSetType(tao->linesearch,morethuente_type);
403:   TaoLineSearchUseTaoRoutines(tao->linesearch,tao);

405:   KSPCreate(((PetscObject)tao)->comm, &tao->ksp);
406:   KSPSetType(tao->ksp,KSPSTCG);
407:   return(0);
408: }
409: EXTERN_C_END