Actual source code: gpcg.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1:  #include <petscksp.h>
  2:  #include <../src/tao/quadratic/impls/gpcg/gpcg.h>


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

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

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

 30: /*------------------------------------------------------------*/
 31: static PetscErrorCode TaoSetFromOptions_GPCG(PetscOptionItems *PetscOptionsObject,Tao tao)
 32: {
 33:   TAO_GPCG       *gpcg = (TAO_GPCG *)tao->data;
 35:   PetscBool      flg;

 38:   PetscOptionsHead(PetscOptionsObject,"Gradient Projection, Conjugate Gradient method for bound constrained optimization");
 39:   ierr=PetscOptionsInt("-tao_gpcg_maxpgits","maximum number of gradient projections per GPCG iterate",NULL,gpcg->maxgpits,&gpcg->maxgpits,&flg);
 40:   PetscOptionsTail();
 41:   KSPSetFromOptions(tao->ksp);
 42:   TaoLineSearchSetFromOptions(tao->linesearch);
 43:   return(0);
 44: }

 46: /*------------------------------------------------------------*/
 47: static PetscErrorCode TaoView_GPCG(Tao tao, PetscViewer viewer)
 48: {
 49:   TAO_GPCG       *gpcg = (TAO_GPCG *)tao->data;
 50:   PetscBool      isascii;

 54:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
 55:   if (isascii) {
 56:     PetscViewerASCIIPrintf(viewer,"Total PG its: %D,",gpcg->total_gp_its);
 57:     PetscViewerASCIIPrintf(viewer,"PG tolerance: %g \n",(double)gpcg->pg_ftol);
 58:   }
 59:   TaoLineSearchView(tao->linesearch,viewer);
 60:   return(0);
 61: }

 63: /* GPCGObjectiveAndGradient()
 64:    Compute f=0.5 * x'Hx + b'x + c
 65:            g=Hx + b
 66: */
 67: static PetscErrorCode GPCGObjectiveAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void*tptr)
 68: {
 69:   Tao            tao = (Tao)tptr;
 70:   TAO_GPCG       *gpcg = (TAO_GPCG*)tao->data;
 72:   PetscReal      f1,f2;

 75:   MatMult(tao->hessian,X,G);
 76:   VecDot(G,X,&f1);
 77:   VecDot(gpcg->B,X,&f2);
 78:   VecAXPY(G,1.0,gpcg->B);
 79:   *f=f1/2.0 + f2 + gpcg->c;
 80:   return(0);
 81: }

 83: /* ---------------------------------------------------------- */
 84: static PetscErrorCode TaoSetup_GPCG(Tao tao)
 85: {
 87:   TAO_GPCG       *gpcg = (TAO_GPCG *)tao->data;

 90:   /* Allocate some arrays */
 91:   if (!tao->gradient) {
 92:       VecDuplicate(tao->solution, &tao->gradient);
 93:   }
 94:   if (!tao->stepdirection) {
 95:       VecDuplicate(tao->solution, &tao->stepdirection);
 96:   }
 97:   if (!tao->XL) {
 98:       VecDuplicate(tao->solution,&tao->XL);
 99:       VecSet(tao->XL,PETSC_NINFINITY);
100:   }
101:   if (!tao->XU) {
102:       VecDuplicate(tao->solution,&tao->XU);
103:       VecSet(tao->XU,PETSC_INFINITY);
104:   }

106:   ierr=VecDuplicate(tao->solution,&gpcg->B);
107:   ierr=VecDuplicate(tao->solution,&gpcg->Work);
108:   ierr=VecDuplicate(tao->solution,&gpcg->X_New);
109:   ierr=VecDuplicate(tao->solution,&gpcg->G_New);
110:   ierr=VecDuplicate(tao->solution,&gpcg->DXFree);
111:   ierr=VecDuplicate(tao->solution,&gpcg->R);
112:   ierr=VecDuplicate(tao->solution,&gpcg->PG);
113:   /*
114:     if (gpcg->ksp_type == GPCG_KSP_NASH) {
115:         KSPSetType(tao->ksp,KSPCGNASH);
116:       } else if (gpcg->ksp_type == GPCG_KSP_STCG) {
117:         KSPSetType(tao->ksp,KSPCGSTCG);
118:       } else {
119:         KSPSetType(tao->ksp,KSPCGGLTR);
120:       }
121:       if (tao->ksp->ops->setfromoptions) {
122:         (*tao->ksp->ops->setfromoptions)(tao->ksp);
123:       }

125:     }
126:   */
127:   return(0);
128: }

130: static PetscErrorCode TaoSolve_GPCG(Tao tao)
131: {
132:   TAO_GPCG                     *gpcg = (TAO_GPCG *)tao->data;
133:   PetscErrorCode               ierr;
134:   PetscInt                     its;
135:   PetscReal                    actred,f,f_new,gnorm,gdx,stepsize,xtb;
136:   PetscReal                    xtHx;
137:   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;


141:   TaoComputeVariableBounds(tao);
142:   VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);
143:   TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);

145:   /* Using f = .5*x'Hx + x'b + c and g=Hx + b,  compute b,c */
146:   TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);
147:   TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);
148:   VecCopy(tao->gradient, gpcg->B);
149:   MatMult(tao->hessian,tao->solution,gpcg->Work);
150:   VecDot(gpcg->Work, tao->solution, &xtHx);
151:   VecAXPY(gpcg->B,-1.0,gpcg->Work);
152:   VecDot(gpcg->B,tao->solution,&xtb);
153:   gpcg->c=f-xtHx/2.0-xtb;
154:   if (gpcg->Free_Local) {
155:       ISDestroy(&gpcg->Free_Local);
156:   }
157:   VecWhichInactive(tao->XL,tao->solution,tao->gradient,tao->XU,PETSC_TRUE,&gpcg->Free_Local);

159:   /* Project the gradient and calculate the norm */
160:   VecCopy(tao->gradient,gpcg->G_New);
161:   VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,gpcg->PG);
162:   VecNorm(gpcg->PG,NORM_2,&gpcg->gnorm);
163:   tao->step=1.0;
164:   gpcg->f = f;

166:     /* Check Stopping Condition      */
167:   tao->reason = TAO_CONTINUE_ITERATING;
168:   TaoLogConvergenceHistory(tao,f,gpcg->gnorm,0.0,tao->ksp_its);
169:   TaoMonitor(tao,tao->niter,f,gpcg->gnorm,0.0,tao->step);
170:   (*tao->ops->convergencetest)(tao,tao->cnvP);

172:   while (tao->reason == TAO_CONTINUE_ITERATING){
173:     tao->ksp_its=0;

175:     GPCGGradProjections(tao);
176:     ISGetSize(gpcg->Free_Local,&gpcg->n_free);

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

180:     KSPReset(tao->ksp);

182:     if (gpcg->n_free > 0){
183:       /* Create a reduced linear system */
184:       VecDestroy(&gpcg->R);
185:       VecDestroy(&gpcg->DXFree);
186:       TaoVecGetSubVec(tao->gradient,gpcg->Free_Local, tao->subset_type, 0.0, &gpcg->R);
187:       VecScale(gpcg->R, -1.0);
188:       TaoVecGetSubVec(tao->stepdirection,gpcg->Free_Local,tao->subset_type, 0.0, &gpcg->DXFree);
189:       VecSet(gpcg->DXFree,0.0);

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

193:       if (tao->hessian_pre == tao->hessian) {
194:         MatDestroy(&gpcg->Hsub_pre);
195:         PetscObjectReference((PetscObject)gpcg->Hsub);
196:         gpcg->Hsub_pre = gpcg->Hsub;
197:       }  else {
198:         TaoMatGetSubMat(tao->hessian, gpcg->Free_Local, gpcg->Work, tao->subset_type, &gpcg->Hsub_pre);
199:       }

201:       KSPReset(tao->ksp);
202:       KSPSetOperators(tao->ksp,gpcg->Hsub,gpcg->Hsub_pre);

204:       KSPSolve(tao->ksp,gpcg->R,gpcg->DXFree);
205:       KSPGetIterationNumber(tao->ksp,&its);
206:       tao->ksp_its+=its;
207:       tao->ksp_tot_its+=its;
208:       VecSet(tao->stepdirection,0.0);
209:       VecISAXPY(tao->stepdirection,gpcg->Free_Local,1.0,gpcg->DXFree);

211:       VecDot(tao->stepdirection,tao->gradient,&gdx);
212:       TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
213:       f_new=f;
214:       TaoLineSearchApply(tao->linesearch,tao->solution,&f_new,tao->gradient,tao->stepdirection,&stepsize,&ls_status);

216:       actred = f_new - f;

218:       /* Evaluate the function and gradient at the new point */
219:       VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU, gpcg->PG);
220:       VecNorm(gpcg->PG, NORM_2, &gnorm);
221:       f=f_new;
222:       ISDestroy(&gpcg->Free_Local);
223:       VecWhichInactive(tao->XL,tao->solution,tao->gradient,tao->XU,PETSC_TRUE,&gpcg->Free_Local);
224:     } else {
225:       actred = 0; gpcg->step=1.0;
226:       /* if there were no free variables, no cg method */
227:     }

229:     tao->niter++;
230:     gpcg->f=f;gpcg->gnorm=gnorm; gpcg->actred=actred;
231:     TaoLogConvergenceHistory(tao,f,gpcg->gnorm,0.0,tao->ksp_its);
232:     TaoMonitor(tao,tao->niter,f,gpcg->gnorm,0.0,tao->step);
233:     (*tao->ops->convergencetest)(tao,tao->cnvP);
234:     if (tao->reason != TAO_CONTINUE_ITERATING) break;
235:   }  /* END MAIN LOOP  */

237:   return(0);
238: }

240: static PetscErrorCode GPCGGradProjections(Tao tao)
241: {
242:   PetscErrorCode                 ierr;
243:   TAO_GPCG                       *gpcg = (TAO_GPCG *)tao->data;
244:   PetscInt                       i;
245:   PetscReal                      actred=-1.0,actred_max=0.0, gAg,gtg=gpcg->gnorm,alpha;
246:   PetscReal                      f_new,gdx,stepsize;
247:   Vec                            DX=tao->stepdirection,XL=tao->XL,XU=tao->XU,Work=gpcg->Work;
248:   Vec                            X=tao->solution,G=tao->gradient;
249:   TaoLineSearchConvergedReason lsflag=TAOLINESEARCH_CONTINUE_ITERATING;

251:   /*
252:      The free, active, and binding variables should be already identified
253:   */
255:   for (i=0;i<gpcg->maxgpits;i++){
256:     if ( -actred <= (gpcg->pg_ftol)*actred_max) break;
257:     VecBoundGradientProjection(G,X,XL,XU,DX);
258:     VecScale(DX,-1.0);
259:     VecDot(DX,G,&gdx);

261:     MatMult(tao->hessian,DX,Work);
262:     VecDot(DX,Work,&gAg);

264:     gpcg->gp_iterates++;
265:     gpcg->total_gp_its++;

267:     gtg=-gdx;
268:     if (PetscAbsReal(gAg) == 0.0) {
269:       alpha = 1.0;
270:     } else {
271:       alpha = PetscAbsReal(gtg/gAg);
272:     }
273:     TaoLineSearchSetInitialStepLength(tao->linesearch,alpha);
274:     f_new=gpcg->f;
275:     TaoLineSearchApply(tao->linesearch,X,&f_new,G,DX,&stepsize,&lsflag);

277:     /* Update the iterate */
278:     actred = f_new - gpcg->f;
279:     actred_max = PetscMax(actred_max,-(f_new - gpcg->f));
280:     gpcg->f = f_new;
281:     ISDestroy(&gpcg->Free_Local);
282:     VecWhichInactive(XL,X,tao->gradient,XU,PETSC_TRUE,&gpcg->Free_Local);
283:   }

285:   gpcg->gnorm=gtg;
286:   return(0);
287: } /* End gradient projections */

289: static PetscErrorCode TaoComputeDual_GPCG(Tao tao, Vec DXL, Vec DXU)
290: {
291:   TAO_GPCG       *gpcg = (TAO_GPCG *)tao->data;

295:   VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, gpcg->Work);
296:   VecCopy(gpcg->Work, DXL);
297:   VecAXPY(DXL,-1.0,tao->gradient);
298:   VecSet(DXU,0.0);
299:   VecPointwiseMax(DXL,DXL,DXU);

301:   VecCopy(tao->gradient,DXU);
302:   VecAXPY(DXU,-1.0,gpcg->Work);
303:   VecSet(gpcg->Work,0.0);
304:   VecPointwiseMin(DXU,gpcg->Work,DXU);
305:   return(0);
306: }

308: /*------------------------------------------------------------*/
309: /*MC
310:   TAOGPCG - gradient projected conjugate gradient algorithm is an active-set
311:         conjugate-gradient based method for bound-constrained minimization

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

317:   Level: beginner
318: M*/
319: PETSC_EXTERN PetscErrorCode TaoCreate_GPCG(Tao tao)
320: {
321:   TAO_GPCG       *gpcg;

325:   tao->ops->setup = TaoSetup_GPCG;
326:   tao->ops->solve = TaoSolve_GPCG;
327:   tao->ops->view  = TaoView_GPCG;
328:   tao->ops->setfromoptions = TaoSetFromOptions_GPCG;
329:   tao->ops->destroy = TaoDestroy_GPCG;
330:   tao->ops->computedual = TaoComputeDual_GPCG;

332:   PetscNewLog(tao,&gpcg);
333:   tao->data = (void*)gpcg;

335:   /* Override default settings (unless already changed) */
336:   if (!tao->max_it_changed) tao->max_it=500;
337:   if (!tao->max_funcs_changed) tao->max_funcs = 100000;
338: #if defined(PETSC_USE_REAL_SINGLE)
339:   if (!tao->gatol_changed) tao->gatol=1e-6;
340:   if (!tao->grtol_changed) tao->grtol=1e-6;
341: #else
342:   if (!tao->gatol_changed) tao->gatol=1e-12;
343:   if (!tao->grtol_changed) tao->grtol=1e-12;
344: #endif

346:   /* Initialize pointers and variables */
347:   gpcg->n=0;
348:   gpcg->maxgpits = 8;
349:   gpcg->pg_ftol = 0.1;

351:   gpcg->gp_iterates=0; /* Cumulative number */
352:   gpcg->total_gp_its = 0;

354:   /* Initialize pointers and variables */
355:   gpcg->n_bind=0;
356:   gpcg->n_free = 0;
357:   gpcg->n_upper=0;
358:   gpcg->n_lower=0;
359:   gpcg->subset_type = TAO_SUBSET_MASK;
360:   gpcg->Hsub=NULL;
361:   gpcg->Hsub_pre=NULL;

363:   KSPCreate(((PetscObject)tao)->comm, &tao->ksp);
364:   PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);
365:   KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);
366:   KSPSetType(tao->ksp,KSPCGNASH);

368:   TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);
369:   PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);
370:   TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHGPCG);
371:   TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch, GPCGObjectiveAndGradient, tao);
372:   TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);
373:   return(0);
374: }