Actual source code: gpcg.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: #include <petsc-private/kspimpl.h>
  2: #include <../src/tao/bound/impls/gpcg/gpcg.h>        /*I "gpcg.h" I*/


  5: static PetscErrorCode GPCGGradProjections(Tao tao);
  6: static PetscErrorCode GPCGObjectiveAndGradient(TaoLineSearch,Vec,PetscReal*,Vec,void*);

  8: /*------------------------------------------------------------*/
 11: static PetscErrorCode TaoDestroy_GPCG(Tao tao)
 12: {
 13:   TAO_GPCG       *gpcg = (TAO_GPCG *)tao->data;

 16:   /* Free allocated memory in GPCG structure */
 18:   VecDestroy(&gpcg->B);
 19:   VecDestroy(&gpcg->Work);
 20:   VecDestroy(&gpcg->X_New);
 21:   VecDestroy(&gpcg->G_New);
 22:   VecDestroy(&gpcg->DXFree);
 23:   VecDestroy(&gpcg->R);
 24:   VecDestroy(&gpcg->PG);
 25:   MatDestroy(&gpcg->Hsub);
 26:   MatDestroy(&gpcg->Hsub_pre);
 27:   ISDestroy(&gpcg->Free_Local);
 28:   PetscFree(tao->data);
 29:   tao->data = NULL;
 30:   return(0);
 31: }

 33: /*------------------------------------------------------------*/
 36: static PetscErrorCode TaoSetFromOptions_GPCG(Tao tao)
 37: {
 38:   TAO_GPCG       *gpcg = (TAO_GPCG *)tao->data;
 40:   PetscBool      flg;

 43:   PetscOptionsHead("Gradient Projection, Conjugate Gradient method for bound constrained optimization");
 44:   ierr=PetscOptionsInt("-tao_gpcg_maxpgits","maximum number of gradient projections per GPCG iterate",0,gpcg->maxgpits,&gpcg->maxgpits,&flg);
 45:   PetscOptionsTail();
 46:   KSPSetFromOptions(tao->ksp);
 47:   TaoLineSearchSetFromOptions(tao->linesearch);
 48:   return(0);
 49: }

 51: /*------------------------------------------------------------*/
 54: static PetscErrorCode TaoView_GPCG(Tao tao, PetscViewer viewer)
 55: {
 56:   TAO_GPCG       *gpcg = (TAO_GPCG *)tao->data;
 57:   PetscBool      isascii;

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

 72: /* GPCGObjectiveAndGradient()
 73:    Compute f=0.5 * x'Hx + b'x + c
 74:            g=Hx + b
 75: */
 78: static PetscErrorCode GPCGObjectiveAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void*tptr)
 79: {
 80:   Tao            tao = (Tao)tptr;
 81:   TAO_GPCG       *gpcg = (TAO_GPCG*)tao->data;
 83:   PetscReal      f1,f2;

 86:   MatMult(tao->hessian,X,G);
 87:   VecDot(G,X,&f1);
 88:   VecDot(gpcg->B,X,&f2);
 89:   VecAXPY(G,1.0,gpcg->B);
 90:   *f=f1/2.0 + f2 + gpcg->c;
 91:   return(0);
 92: }

 94: /* ---------------------------------------------------------- */
 97: static PetscErrorCode TaoSetup_GPCG(Tao tao)
 98: {
100:   TAO_GPCG       *gpcg = (TAO_GPCG *)tao->data;

103:   /* Allocate some arrays */
104:   if (!tao->gradient) {
105:       VecDuplicate(tao->solution, &tao->gradient);
106:   }
107:   if (!tao->stepdirection) {
108:       VecDuplicate(tao->solution, &tao->stepdirection);
109:   }
110:   if (!tao->XL) {
111:       VecDuplicate(tao->solution,&tao->XL);
112:       VecSet(tao->XL,PETSC_NINFINITY);
113:   }
114:   if (!tao->XU) {
115:       VecDuplicate(tao->solution,&tao->XU);
116:       VecSet(tao->XU,PETSC_INFINITY);
117:   }

119:   ierr=VecDuplicate(tao->solution,&gpcg->B);
120:   ierr=VecDuplicate(tao->solution,&gpcg->Work);
121:   ierr=VecDuplicate(tao->solution,&gpcg->X_New);
122:   ierr=VecDuplicate(tao->solution,&gpcg->G_New);
123:   ierr=VecDuplicate(tao->solution,&gpcg->DXFree);
124:   ierr=VecDuplicate(tao->solution,&gpcg->R);
125:   ierr=VecDuplicate(tao->solution,&gpcg->PG);
126:   /*
127:     if (gpcg->ksp_type == GPCG_KSP_NASH) {
128:         KSPSetType(tao->ksp,KSPNASH);
129:       } else if (gpcg->ksp_type == GPCG_KSP_STCG) {
130:         KSPSetType(tao->ksp,KSPSTCG);
131:       } else {
132:         KSPSetType(tao->ksp,KSPGLTR);
133:       }
134:       if (tao->ksp->ops->setfromoptions) {
135:         (*tao->ksp->ops->setfromoptions)(tao->ksp);
136:       }

138:     }
139:   */
140:   return(0);
141: }

145: static PetscErrorCode TaoSolve_GPCG(Tao tao)
146: {
147:   TAO_GPCG                     *gpcg = (TAO_GPCG *)tao->data;
148:   PetscErrorCode               ierr;
149:   PetscInt                     iter=0,its;
150:   PetscReal                    actred,f,f_new,gnorm,gdx,stepsize,xtb;
151:   PetscReal                    xtHx;
152:   TaoConvergedReason           reason = TAO_CONTINUE_ITERATING;
153:   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;

156:   gpcg->Hsub=NULL;
157:   gpcg->Hsub_pre=NULL;

159:   TaoComputeVariableBounds(tao);
160:   VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);
161:   TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);

163:   /* Using f = .5*x'Hx + x'b + c and g=Hx + b,  compute b,c */
164:   TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);
165:   TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);
166:   VecCopy(tao->gradient, gpcg->B);
167:   MatMult(tao->hessian,tao->solution,gpcg->Work);
168:   VecDot(gpcg->Work, tao->solution, &xtHx);
169:   VecAXPY(gpcg->B,-1.0,gpcg->Work);
170:   VecDot(gpcg->B,tao->solution,&xtb);
171:   gpcg->c=f-xtHx/2.0-xtb;
172:   if (gpcg->Free_Local) {
173:       ISDestroy(&gpcg->Free_Local);
174:   }
175:   VecWhichBetween(tao->XL,tao->solution,tao->XU,&gpcg->Free_Local);

177:   /* Project the gradient and calculate the norm */
178:   VecCopy(tao->gradient,gpcg->G_New);
179:   VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,gpcg->PG);
180:   VecNorm(gpcg->PG,NORM_2,&gpcg->gnorm);
181:   tao->step=1.0;
182:   gpcg->f = f;

184:     /* Check Stopping Condition      */
185:   ierr=TaoMonitor(tao,iter,f,gpcg->gnorm,0.0,tao->step,&reason);

187:   while (reason == TAO_CONTINUE_ITERATING){

189:     GPCGGradProjections(tao);
190:     ISGetSize(gpcg->Free_Local,&gpcg->n_free);

192:     f=gpcg->f; gnorm=gpcg->gnorm;

194:     KSPReset(tao->ksp);

196:     if (gpcg->n_free > 0){
197:       /* Create a reduced linear system */
198:       VecDestroy(&gpcg->R);
199:       VecDestroy(&gpcg->DXFree);
200:       TaoVecGetSubVec(tao->gradient,gpcg->Free_Local, tao->subset_type, 0.0, &gpcg->R);
201:       VecScale(gpcg->R, -1.0);
202:       TaoVecGetSubVec(tao->stepdirection,gpcg->Free_Local,tao->subset_type, 0.0, &gpcg->DXFree);
203:       VecSet(gpcg->DXFree,0.0);

205:       TaoMatGetSubMat(tao->hessian, gpcg->Free_Local, gpcg->Work, tao->subset_type, &gpcg->Hsub);

207:       if (tao->hessian_pre == tao->hessian) {
208:         MatDestroy(&gpcg->Hsub_pre);
209:         PetscObjectReference((PetscObject)gpcg->Hsub);
210:         gpcg->Hsub_pre = gpcg->Hsub;
211:       }  else {
212:         TaoMatGetSubMat(tao->hessian, gpcg->Free_Local, gpcg->Work, tao->subset_type, &gpcg->Hsub_pre);
213:       }

215:       KSPReset(tao->ksp);
216:       KSPSetOperators(tao->ksp,gpcg->Hsub,gpcg->Hsub_pre);

218:       KSPSolve(tao->ksp,gpcg->R,gpcg->DXFree);
219:       KSPGetIterationNumber(tao->ksp,&its);
220:       tao->ksp_its+=its;

222:       VecSet(tao->stepdirection,0.0);
223:       VecISAXPY(tao->stepdirection,gpcg->Free_Local,1.0,gpcg->DXFree);

225:       VecDot(tao->stepdirection,tao->gradient,&gdx);
226:       TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
227:       f_new=f;
228:       TaoLineSearchApply(tao->linesearch,tao->solution,&f_new,tao->gradient,tao->stepdirection,&stepsize,&ls_status);

230:       actred = f_new - f;

232:       /* Evaluate the function and gradient at the new point */
233:       VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU, gpcg->PG);
234:       VecNorm(gpcg->PG, NORM_2, &gnorm);
235:       f=f_new;
236:       ISDestroy(&gpcg->Free_Local);
237:       VecWhichBetween(tao->XL,tao->solution,tao->XU,&gpcg->Free_Local);
238:     } else {
239:       actred = 0; gpcg->step=1.0;
240:       /* if there were no free variables, no cg method */
241:     }

243:     iter++;
244:     TaoMonitor(tao,iter,f,gnorm,0.0,gpcg->step,&reason);
245:     gpcg->f=f;gpcg->gnorm=gnorm; gpcg->actred=actred;
246:     if (reason!=TAO_CONTINUE_ITERATING) break;
247:   }  /* END MAIN LOOP  */

249:   return(0);
250: }

254: static PetscErrorCode GPCGGradProjections(Tao tao)
255: {
256:   PetscErrorCode                 ierr;
257:   TAO_GPCG                       *gpcg = (TAO_GPCG *)tao->data;
258:   PetscInt                       i;
259:   PetscReal                      actred=-1.0,actred_max=0.0, gAg,gtg=gpcg->gnorm,alpha;
260:   PetscReal                      f_new,gdx,stepsize;
261:   Vec                            DX=tao->stepdirection,XL=tao->XL,XU=tao->XU,Work=gpcg->Work;
262:   Vec                            X=tao->solution,G=tao->gradient;
263:   TaoLineSearchConvergedReason lsflag=TAOLINESEARCH_CONTINUE_ITERATING;

265:   /*
266:      The free, active, and binding variables should be already identified
267:   */
269:   for (i=0;i<gpcg->maxgpits;i++){
270:     if ( -actred <= (gpcg->pg_ftol)*actred_max) break;
271:     VecBoundGradientProjection(G,X,XL,XU,DX);
272:     VecScale(DX,-1.0);
273:     VecDot(DX,G,&gdx);

275:     MatMult(tao->hessian,DX,Work);
276:     VecDot(DX,Work,&gAg);

278:     gpcg->gp_iterates++;
279:     gpcg->total_gp_its++;

281:     gtg=-gdx;
282:     alpha = PetscAbsReal(gtg/gAg);
283:     TaoLineSearchSetInitialStepLength(tao->linesearch,alpha);
284:     f_new=gpcg->f;
285:     TaoLineSearchApply(tao->linesearch,X,&f_new,G,DX,&stepsize,&lsflag);

287:     /* Update the iterate */
288:     actred = f_new - gpcg->f;
289:     actred_max = PetscMax(actred_max,-(f_new - gpcg->f));
290:     gpcg->f = f_new;
291:     ISDestroy(&gpcg->Free_Local);
292:     VecWhichBetween(XL,X,XU,&gpcg->Free_Local);
293:   }

295:   gpcg->gnorm=gtg;
296:   return(0);
297: } /* End gradient projections */

301: static PetscErrorCode TaoComputeDual_GPCG(Tao tao, Vec DXL, Vec DXU)
302: {
303:   TAO_GPCG       *gpcg = (TAO_GPCG *)tao->data;

307:   VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, gpcg->Work);
308:   VecCopy(gpcg->Work, DXL);
309:   VecAXPY(DXL,-1.0,tao->gradient);
310:   VecSet(DXU,0.0);
311:   VecPointwiseMax(DXL,DXL,DXU);

313:   VecCopy(tao->gradient,DXU);
314:   VecAXPY(DXU,-1.0,gpcg->Work);
315:   VecSet(gpcg->Work,0.0);
316:   VecPointwiseMin(DXU,gpcg->Work,DXU);
317:   return(0);
318: }

320: /*------------------------------------------------------------*/
321: /*MC
322:   TAOGPCG - gradient projected conjugate gradient algorithm is an active-set
323:         conjugate-gradient based method for bound-constrained minimization

325:   Options Database Keys:
326: + -tao_gpcg_maxpgits - maximum number of gradient projections for GPCG iterate
327: - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets

329:   Level: beginner
330: M*/
331: EXTERN_C_BEGIN
334: PetscErrorCode TaoCreate_GPCG(Tao tao)
335: {
336:   TAO_GPCG       *gpcg;

340:   tao->ops->setup = TaoSetup_GPCG;
341:   tao->ops->solve = TaoSolve_GPCG;
342:   tao->ops->view  = TaoView_GPCG;
343:   tao->ops->setfromoptions = TaoSetFromOptions_GPCG;
344:   tao->ops->destroy = TaoDestroy_GPCG;
345:   tao->ops->computedual = TaoComputeDual_GPCG;

347:   PetscNewLog(tao,&gpcg);
348:   tao->data = (void*)gpcg;

350:   tao->max_it = 500;
351:   tao->max_funcs = 100000;

353: #if defined(PETSC_USE_REAL_SINGLE)
354:   tao->fatol = 1e-6;
355:   tao->frtol = 1e-6;
356: #else
357:   tao->fatol = 1e-12;
358:   tao->frtol = 1e-12;
359: #endif

361:   /* Initialize pointers and variables */
362:   gpcg->n=0;
363:   gpcg->maxgpits = 8;
364:   gpcg->pg_ftol = 0.1;

366:   gpcg->gp_iterates=0; /* Cumulative number */
367:   gpcg->total_gp_its = 0;

369:   /* Initialize pointers and variables */
370:   gpcg->n_bind=0;
371:   gpcg->n_free = 0;
372:   gpcg->n_upper=0;
373:   gpcg->n_lower=0;
374:   gpcg->subset_type = TAO_SUBSET_MASK;
375:   /* gpcg->ksp_type = GPCG_KSP_STCG; */

377:   KSPCreate(((PetscObject)tao)->comm, &tao->ksp);
378:   KSPSetType(tao->ksp,KSPNASH);

380:   TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);
381:   TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHGPCG);
382:   TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch, GPCGObjectiveAndGradient, tao);
383:   return(0);
384: }
385: EXTERN_C_END