Actual source code: tron.c

petsc-3.7.5 2017-01-01
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(PetscOptionItems *PetscOptionsObject,Tao tao)
 36: {
 37:   TAO_TRON       *tron = (TAO_TRON *)tao->data;
 39:   PetscBool      flg;

 42:   PetscOptionsHead(PetscOptionsObject,"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                     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, tao->niter, tron->f, tron->gnorm, 0.0, tron->stepsize, &reason);
136:   while (reason==TAO_CONTINUE_ITERATING){
137:     tao->ksp_its=0;
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:       PetscInfo(tao,"No free variables in tron iteration.\n");
148:       VecNorm(tao->gradient,NORM_2,&tron->gnorm);
149:       TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, delta, &reason);
150:       if (!reason) {
151:         reason = TAO_CONVERGED_STEPTOL;
152:         TaoSetConvergedReason(tao,reason);
153:       }
154:       break;
155:     }
156:     /* use free_local to mask/submat gradient, hessian, stepdirection */
157:     TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->R);
158:     TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->DXFree);
159:     VecSet(tron->DXFree,0.0);
160:     VecScale(tron->R, -1.0);
161:     TaoMatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub);
162:     if (tao->hessian == tao->hessian_pre) {
163:       MatDestroy(&tron->Hpre_sub);
164:       PetscObjectReference((PetscObject)(tron->H_sub));
165:       tron->Hpre_sub = tron->H_sub;
166:     } else {
167:       TaoMatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type,&tron->Hpre_sub);
168:     }
169:     KSPReset(tao->ksp);
170:     KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub);
171:     while (1) {

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

176:       KSPSolve(tao->ksp, tron->R, tron->DXFree);
177:       KSPGetIterationNumber(tao->ksp,&its);
178:       tao->ksp_its+=its;
179:       tao->ksp_tot_its+=its;
180:       VecSet(tao->stepdirection,0.0);

182:       /* Add dxfree matrix to compute step direction vector */
183:       VecISAXPY(tao->stepdirection,tron->Free_Local,1.0,tron->DXFree);

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

188:       VecCopy(tao->solution, tron->X_New);
189:       VecCopy(tao->gradient, tron->G_New);

191:       stepsize=1.0;f_new=f;

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

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

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

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

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

243:     tron->f=f; tron->actred=actred; tao->trust=delta;
244:     tao->niter++;
245:     TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, delta, &reason);
246:   }  /* END MAIN LOOP  */

248:   return(0);
249: }


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

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

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

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

277:     tron->gp_iterates++; tron->total_gp_its++;
278:     f_new=tron->f;

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


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

298:   return(0);
299: }

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

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

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

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

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

327: /*------------------------------------------------------------*/
328: /*MC
329:   TAOTRON - The TRON algorithm is an active-set Newton trust region method
330:   for bound-constrained minimization.

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

336:   Level: beginner
337: M*/
340: PETSC_EXTERN PetscErrorCode TaoCreate_TRON(Tao tao)
341: {
342:   TAO_TRON       *tron;
344:   const char     *morethuente_type = TAOLINESEARCHMT;

347:   tao->ops->setup = TaoSetup_TRON;
348:   tao->ops->solve = TaoSolve_TRON;
349:   tao->ops->view = TaoView_TRON;
350:   tao->ops->setfromoptions = TaoSetFromOptions_TRON;
351:   tao->ops->destroy = TaoDestroy_TRON;
352:   tao->ops->computedual = TaoComputeDual_TRON;

354:   PetscNewLog(tao,&tron);
355:   tao->data = (void*)tron;

357:   /* Override default settings (unless already changed) */
358:   if (!tao->max_it_changed) tao->max_it = 50;

360: #if defined(PETSC_USE_REAL_SINGLE)
361:   if (!tao->steptol_changed) tao->steptol = 1.0e-6;
362: #else
363:   if (!tao->steptol_changed) tao->steptol = 1.0e-12;
364: #endif

366:   if (!tao->trust0_changed) tao->trust0 = 1.0;

368:   /* Initialize pointers and variables */
369:   tron->n            = 0;
370:   tron->maxgpits     = 3;
371:   tron->pg_ftol      = 0.001;

373:   tron->eta1         = 1.0e-4;
374:   tron->eta2         = 0.25;
375:   tron->eta3         = 0.50;
376:   tron->eta4         = 0.90;

378:   tron->sigma1       = 0.5;
379:   tron->sigma2       = 2.0;
380:   tron->sigma3       = 4.0;

382:   tron->gp_iterates  = 0; /* Cumulative number */
383:   tron->total_gp_its = 0;
384:   tron->n_free       = 0;

386:   tron->DXFree=NULL;
387:   tron->R=NULL;
388:   tron->X_New=NULL;
389:   tron->G_New=NULL;
390:   tron->Work=NULL;
391:   tron->Free_Local=NULL;
392:   tron->H_sub=NULL;
393:   tron->Hpre_sub=NULL;
394:   tao->subset_type = TAO_SUBSET_SUBVEC;

396:   TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);
397:   TaoLineSearchSetType(tao->linesearch,morethuente_type);
398:   TaoLineSearchUseTaoRoutines(tao->linesearch,tao);
399:   TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);

401:   KSPCreate(((PetscObject)tao)->comm, &tao->ksp);
402:   KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);
403:   KSPSetType(tao->ksp,KSPSTCG);
404:   return(0);
405: }