Actual source code: gpcg.c

  1: /*$Id$*/

  3: #include "gpcg.h"        /*I "tao_solver.h" I*/

  5: char gpcgname[]="tao_gpcg";
  6: TaoMethod gpcgtypename = gpcgname;

  8: static int TaoGradProjections(TAO_SOLVER, TAO_GPCG *);
  9: static int GPCGCheckOptimalFace(TaoVec*, TaoVec*, TaoVec*, TaoVec*, TaoVec*, 
 10:                                 TaoIndexSet*, TaoIndexSet*, TaoTruth *optimal);

 12: /*------------------------------------------------------------*/
 15: static int TaoSetDown_GPCG(TAO_SOLVER tao, void*solver)
 16: {
 17:   TAO_GPCG *gpcg = (TAO_GPCG *)solver;
 18:   int      info;
 19:   /* Free allocated memory in GPCG structure */
 20:   TaoFunctionBegin;
 21:   
 22:   info = TaoVecDestroy(gpcg->X_New);CHKERRQ(info);
 23:   info = TaoVecDestroy(gpcg->G_New);CHKERRQ(info);
 24:   info = TaoVecDestroy(gpcg->DX);CHKERRQ(info);gpcg->DX=0;
 25:   info = TaoVecDestroy(gpcg->Work);CHKERRQ(info);
 26:   info = TaoVecDestroy(gpcg->DXFree);CHKERRQ(info);
 27:   info = TaoVecDestroy(gpcg->R);CHKERRQ(info);
 28:   info = TaoVecDestroy(gpcg->B);CHKERRQ(info);
 29:   info = TaoVecDestroy(gpcg->G);CHKERRQ(info);
 30:   info = TaoVecDestroy(gpcg->PG);CHKERRQ(info);
 31:   info = TaoVecDestroy(gpcg->XL);CHKERRQ(info);
 32:   info = TaoVecDestroy(gpcg->XU);CHKERRQ(info);
 33:   
 34:   info = TaoIndexSetDestroy(gpcg->Free_Local);CHKERRQ(info);
 35:   info = TaoIndexSetDestroy(gpcg->TT);CHKERRQ(info);
 36:   info = TaoMatDestroy(gpcg->Hsub);CHKERRQ(info);
 37:   
 38:   TaoFunctionReturn(0);
 39: }

 41: /*------------------------------------------------------------*/
 44: static int TaoSetFromOptions_GPCG(TAO_SOLVER tao, void*solver)
 45: {
 46:   TAO_GPCG *gpcg = (TAO_GPCG *)solver;
 47:   int      info,ival;
 48:   TaoTruth flg;

 50:   TaoFunctionBegin;
 51:   info = TaoOptionsHead("Gradient Projection, Conjugate Gradient method for bound constrained optimization");CHKERRQ(info);

 53:   info=TaoOptionInt("-gpcg_maxpgits","maximum number of gradient projections per GPCG iterate",0,gpcg->maxgpits,&gpcg->maxgpits,&flg);
 54:   CHKERRQ(info);

 56:   info = TaoOptionInt("-redistribute","Redistribute Free variables (> 1 processors, only)","TaoPetscISType",1,&ival,&flg); CHKERRQ(info);

 58:   info = TaoOptionName("-submatrixfree","Mask full matrix instead of extract submatrices","TaoPetscISType",&flg); CHKERRQ(info);


 61:   info = TaoOptionsTail();CHKERRQ(info);
 62:   info=TaoLineSearchSetFromOptions(tao);CHKERRQ(info);

 64:   TaoFunctionReturn(0);
 65: }

 67: /*------------------------------------------------------------*/
 70: static int TaoView_GPCG(TAO_SOLVER tao,void*solver)
 71: {
 72:   TAO_GPCG *gpcg = (TAO_GPCG *)solver;
 73:   int      info;

 75:   TaoFunctionBegin;

 77:   info = TaoPrintInt(tao," Total PG its: %d,",gpcg->total_gp_its);CHKERRQ(info);
 78:   info = TaoPrintDouble(tao," PG tolerance: %4.3f \n",gpcg->pg_ftol);CHKERRQ(info);
 79:   info = TaoLineSearchView(tao);CHKERRQ(info);

 81:   TaoFunctionReturn(0);
 82: }

 84: /* ---------------------------------------------------------- */
 87: int TaoGPCGComputeFunctionGradient(TAO_SOLVER tao, TaoVec *XX, double *ff, TaoVec *GG){
 88:   TAO_GPCG *gpcg;
 89:   TaoMat *HH;
 90:   int info;
 91:   double f1,f2;

 93:   TaoFunctionBegin;
 94:   info = TaoGetSolverContext(tao,gpcgtypename,(void**)&gpcg); CHKERRQ(info);
 95:   if (gpcg==0){TaoFunctionReturn(0);}
 96:   info = TaoGetHessian(tao,&HH);CHKERRQ(info);
 97:   info = HH->Multiply(XX,GG);CHKERRQ(info);
 98:   info = XX->Dot(GG,&f1);CHKERRQ(info);
 99:   info = XX->Dot(gpcg->B,&f2);CHKERRQ(info);
100:   info = GG->Axpy(1.0,gpcg->B);CHKERRQ(info);
101:   *ff=f1/2.0 + f2 + gpcg->c;
102:   TaoFunctionReturn(0);
103: }

105: /* ---------------------------------------------------------- */
108: static int TaoSetUp_GPCG(TAO_SOLVER tao,void*solver){

110:   int      n,info;
111:   TAO_GPCG *gpcg = (TAO_GPCG *) solver;
112:   TaoVec   *X;
113:   TaoMat   *HH;
114:   TaoIndexSet *TIS;

116:   TaoFunctionBegin;

118:   info = TaoGetSolution(tao,&X);CHKERRQ(info); gpcg->X=X;
119:   info = TaoGetHessian(tao,&HH);CHKERRQ(info); gpcg->H=HH;

121:   /* Allocate some arrays */
122:   info=X->Clone(&gpcg->DX); CHKERRQ(info);
123:   info=X->Clone(&gpcg->B); CHKERRQ(info);
124:   info=X->Clone(&gpcg->Work); CHKERRQ(info);
125:   info=X->Clone(&gpcg->X_New); CHKERRQ(info);
126:   info=X->Clone(&gpcg->G_New); CHKERRQ(info);
127:   info=X->Clone(&gpcg->DXFree); CHKERRQ(info);
128:   info=X->Clone(&gpcg->R); CHKERRQ(info);
129:   info=X->Clone(&gpcg->G); CHKERRQ(info);
130:   info=X->Clone(&gpcg->PG); CHKERRQ(info);
131:   info=X->Clone(&gpcg->XL); CHKERRQ(info);
132:   info=X->Clone(&gpcg->XU); CHKERRQ(info);

134:   info = TaoSetLagrangianGradientVector(tao,gpcg->PG);CHKERRQ(info);
135:   info = TaoSetStepDirectionVector(tao,gpcg->DX);CHKERRQ(info);
136:   info = TaoSetVariableBounds(tao,gpcg->XL,gpcg->XU);CHKERRQ(info);

138:   info = X->GetDimension(&n); CHKERRQ(info);
139:   gpcg->n=n;
140:   info = TaoCreateLinearSolver(tao,HH,300,0); CHKERRQ(info);

142:   info = X->CreateIndexSet(&TIS); CHKERRQ(info);
143:   gpcg->Free_Local = TIS;
144:   info = gpcg->Free_Local->Duplicate(&gpcg->TT); CHKERRQ(info);

146:   info = HH->CreateReducedMatrix(TIS,TIS,&gpcg->Hsub); CHKERRQ(info);

148:   TaoFunctionReturn(0);
149: }

153: static int TaoSolve_GPCG(TAO_SOLVER tao, void *solver)
154: {
155:   TAO_GPCG *gpcg = (TAO_GPCG *)solver ;
156:   int info,lsflag,iter=0;
157:   TaoTruth optimal_face=TAO_FALSE,success;
158:   double actred,f,f_new,f_full,gnorm,gdx,stepsize;
159:   double c;
160:   TaoVec *XU, *XL;
161:   TaoVec *X,  *G=gpcg->G , *B=gpcg->B, *PG=gpcg->PG;
162:   TaoVec *R=gpcg->R, *DXFree=gpcg->DXFree;
163:   TaoVec *G_New=gpcg->G_New;
164:   TaoVec *DX=gpcg->DX, *Work=gpcg->Work;
165:   TaoMat *H, *Hsub=gpcg->Hsub;
166:   TaoIndexSet *Free_Local = gpcg->Free_Local, *TIS=gpcg->TT;
167:   TaoTerminateReason reason;

169:   TaoFunctionBegin;

171:   /* Check if upper bound greater than lower bound. */
172:   info = TaoGetSolution(tao,&X);CHKERRQ(info);
173:   info = TaoGetHessian(tao,&H);CHKERRQ(info);

175:   info = TaoGetVariableBounds(tao,&XL,&XU);CHKERRQ(info);
176:   info = TaoEvaluateVariableBounds(tao,XL,XU); CHKERRQ(info);
177:   info = X->Median(XL,X,XU); CHKERRQ(info);

179:   info = TaoComputeHessian(tao,X,H); CHKERRQ(info);
180:   info = TaoComputeFunctionGradient(tao,X,&f,B);
181:   CHKERRQ(info);

183:   /* Compute quadratic representation */
184:   info = H->Multiply(X,Work); CHKERRQ(info);
185:   info = X->Dot(Work,&c); CHKERRQ(info);
186:   info = B->Axpy(-1.0,Work); CHKERRQ(info);
187:   info = X->Dot(B,&stepsize); CHKERRQ(info);
188:   gpcg->c=f-c/2.0-stepsize;

190:   info = Free_Local->WhichBetween(XL,X,XU); CHKERRQ(info);
191:   
192:   info = TaoGPCGComputeFunctionGradient(tao, X, &gpcg->f , G);
193:   
194:   /* Project the gradient and calculate the norm */
195:   info = G_New->CopyFrom(G);CHKERRQ(info);
196:   info = PG->BoundGradientProjection(G,XL,X,XU);CHKERRQ(info);
197:   info = PG->Norm2(&gpcg->gnorm); CHKERRQ(info);
198:   gpcg->step=1.0;

200:     /* Check Stopping Condition      */
201:   info=TaoMonitor(tao,iter++,gpcg->f,gpcg->gnorm,0,gpcg->step,&reason); CHKERRQ(info);

203:   while (reason == TAO_CONTINUE_ITERATING){

205:     info = TaoGradProjections(tao, gpcg); CHKERRQ(info);

207:     info = Free_Local->WhichBetween(XL,X,XU); CHKERRQ(info);
208:     info = Free_Local->GetSize(&gpcg->n_free); CHKERRQ(info);
209:     f=gpcg->f; gnorm=gpcg->gnorm; 

211:     if (gpcg->n_free > 0){
212:       
213:       /* Create a reduced linear system */
214:       info = R->SetReducedVec(G,Free_Local);CHKERRQ(info);
215:       info = R->Negate(); CHKERRQ(info);
216:       info = DXFree->SetReducedVec(DX,Free_Local);CHKERRQ(info);
217:       info = DXFree->SetToZero(); CHKERRQ(info);

219:       info = Hsub->SetReducedMatrix(H,Free_Local,Free_Local);CHKERRQ(info);

221:       info = TaoPreLinearSolve(tao,Hsub);CHKERRQ(info);

223:       /* Approximately solve the reduced linear system */
224:       info = TaoLinearSolve(tao,Hsub,R,DXFree,&success);CHKERRQ(info);
225:       
226:       info=DX->SetToZero(); CHKERRQ(info);
227:       info=DX->ReducedXPY(DXFree,Free_Local);CHKERRQ(info);
228:       
229:       info = G->Dot(DX,&gdx); CHKERRQ(info);
230:       
231:       stepsize=1.0; f_new=f;
232:       info = TaoLineSearchApply(tao,X,G,DX,Work,
233:                                 &f_new,&f_full,&stepsize,&lsflag);
234:       CHKERRQ(info);
235:       
236:       actred = f_new - f;
237:       
238:       /* Evaluate the function and gradient at the new point */      
239:       info =  PG->BoundGradientProjection(G,XL,X,XU);
240:       CHKERRQ(info);
241:       info = PG->Norm2(&gnorm);  CHKERRQ(info);      
242:       f=f_new;
243:       
244:       info = GPCGCheckOptimalFace(X,XL,XU,PG,Work, Free_Local, TIS,
245:                                   &optimal_face); CHKERRQ(info);
246:       
247:     } else {
248:       
249:       actred = 0; stepsize=1.0;
250:       /* if there were no free variables, no cg method */

252:     }

254:     info = TaoMonitor(tao,iter,f,gnorm,0.0,stepsize,&reason); CHKERRQ(info);
255:     gpcg->f=f;gpcg->gnorm=gnorm; gpcg->actred=actred;
256:     if (reason!=TAO_CONTINUE_ITERATING) break;
257:     iter++;


260:   }  /* END MAIN LOOP  */

262:   TaoFunctionReturn(0);
263: }

267: static int TaoGradProjections(TAO_SOLVER tao, TAO_GPCG *gpcg)
268: {
269:   int i,lsflag=0,info;
270:   TaoTruth optimal_face=TAO_FALSE;
271:   double actred=-1.0,actred_max=0.0, gAg,gtg=gpcg->gnorm,alpha;
272:   double f_new,f_full,gdx;
273:   TaoMat *H;
274:   TaoVec *DX=gpcg->DX,*XL=gpcg->XL,*XU=gpcg->XU,*Work=gpcg->Work;
275:   TaoVec *X=gpcg->X,*G=gpcg->G;
276:   /*
277:      The gradient and function value passed into and out of this
278:      routine should be current and correct.
279:      
280:      The free, active, and binding variables should be already identified
281:   */
282:   
283:   TaoFunctionBegin;
284:   
285:   info = TaoGetSolution(tao,&X);CHKERRQ(info);
286:   info = TaoGetHessian(tao,&H);CHKERRQ(info);
287:   info = TaoGetVariableBounds(tao,&XL,&XU);CHKERRQ(info);

289:   for (i=0;i<gpcg->maxgpits;i++){

291:     if ( -actred <= (gpcg->pg_ftol)*actred_max) break;
292:  
293:     info = DX->BoundGradientProjection(G,XL,X,XU); CHKERRQ(info);
294:     info = DX->Negate(); CHKERRQ(info);
295:     info = DX->Dot(G,&gdx); CHKERRQ(info);

297:     info= H->Multiply(DX,Work); CHKERRQ(info);
298:     info= DX->Dot(Work,&gAg); CHKERRQ(info);
299:  
300:     gpcg->gp_iterates++; gpcg->total_gp_its++;    
301:   
302:     gtg=-gdx;
303:     alpha = TaoAbsDouble(gtg/gAg);
304:     gpcg->stepsize = alpha; f_new=gpcg->f;

306:     info = TaoLineSearchApply(tao,X,G,DX,Work,
307:                               &f_new,&f_full,&gpcg->stepsize,&lsflag);
308:     CHKERRQ(info);

310:     /* Update the iterate */
311:     actred = f_new - gpcg->f;
312:     actred_max = TaoMax(actred_max,-(f_new - gpcg->f));
313:     gpcg->f = f_new;
314:     info = GPCGCheckOptimalFace(X,XL,XU,G,Work,gpcg->Free_Local,gpcg->TT,
315:                                 &optimal_face); CHKERRQ(info);

317:     if ( optimal_face == TAO_TRUE ) break;

319:   }
320:   
321:   gpcg->gnorm=gtg;
322:   TaoFunctionReturn(0);

324: } /* End gradient projections */


329: static int GPCGCheckOptimalFace(TaoVec *X,TaoVec *XL,TaoVec*XU,TaoVec *PG,TaoVec*W,
330:                                 TaoIndexSet*Free_Local, TaoIndexSet*TT,
331:                                 TaoTruth *optimal)
332: {
333:   int info,n_free;
334:   double rr;
335:   TaoTruth same;

337:   TaoFunctionBegin;
338:   *optimal = TAO_FALSE;

340:   /* Also check to see if the active set is the same */

342:   info = TT->WhichBetween(XL,X,XU); CHKERRQ(info);
343:   info = Free_Local->IsSame(TT,&same); CHKERRQ(info);
344:   info = Free_Local->GetSize(&n_free); CHKERRQ(info);
345:   if (same == TAO_FALSE){
346:     info = Free_Local->WhichBetween(XL,X,XU); CHKERRQ(info);
347:     *optimal = TAO_FALSE;
348:     TaoFunctionReturn(0);
349:   } else {
350:     *optimal = TAO_TRUE;
351:   }


354:   info = W->CopyFrom(PG); CHKERRQ(info);
355:   info = W->Negate(); CHKERRQ(info);

357:   info = W->BoundGradientProjection(W,XL,X,XU); CHKERRQ(info);
358:   info = W->Axpy(1.0,PG); CHKERRQ(info);

360:   info = W->Norm2(&rr); CHKERRQ(info);
361:   if (rr>0) *optimal = TAO_FALSE;

363:   *optimal = TAO_FALSE;
364:   /*
365:   info = gpcg->TT->whichNonNegative(W); CHKERRQ(info);
366:   info = gpcg->TT->GetSize(&n); CHKERRQ(info);
367:   if (n==0) *optimal = TAO_TRUE;
368:   */
369:   TaoFunctionReturn(0);
370: }

372: int TaoDefaultMonitor_GPCG(TAO_SOLVER tao,void *dummy)
373: {
374:   TAO_GPCG *gpcg;
375:   double   fct,gnorm;
376:   int      its,nfree,info;

378:   TaoFunctionBegin;
379:   info = TaoGetSolutionStatus(tao,&its,&fct,&gnorm,0,0,0);CHKERRQ(info);
380:   info = TaoGetSolverContext(tao,gpcgtypename,(void**)&gpcg); CHKERRQ(info);
381:   if (gpcg==0){TaoFunctionReturn(0);}
382:   nfree=gpcg->n_free;
383:   info = TaoPrintInt(tao,"iter: %d,",its);CHKERRQ(info);
384:   info = TaoPrintDouble(tao," Fcn value: %g,",fct);CHKERRQ(info);
385:   info = TaoPrintDouble(tao," PGrad. norm: %g, ",gnorm);CHKERRQ(info);
386:   info = TaoPrintInt(tao,"free vars:%d \n",nfree);CHKERRQ(info);
387:   TaoFunctionReturn(0);
388: }

392: int TaoGetDualVariables_GPCG(TAO_SOLVER tao, TaoVec* DXL, TaoVec* DXU, void* solver)
393: {

395:   TAO_GPCG *gpcg = (TAO_GPCG *) solver;
396:   TaoVec  *G=gpcg->G,*GP=gpcg->Work;
397:   TaoVec  *X,*XL,*XU;
398:   int       info;

400:   TaoFunctionBegin;
401:   info = TaoGetVariableBounds(tao,&XL,&XU); CHKERRQ(info);
402:   info = TaoGetSolution(tao,&X); CHKERRQ(info);
403:   info = GP->BoundGradientProjection(G,XL,X,XU); CHKERRQ(info);

405:   info = DXL->Waxpby(-1.0,G,1.0,GP); CHKERRQ(info);
406:   info = DXU->SetToZero(); CHKERRQ(info);
407:   info = DXL->PointwiseMaximum(DXL,DXU); CHKERRQ(info);

409:   info = DXU->Waxpby(-1.0,GP,1.0,G); CHKERRQ(info);
410:   info = GP->SetToZero(); CHKERRQ(info);
411:   info = DXU->PointwiseMinimum(GP,DXU); CHKERRQ(info);

413:   TaoFunctionReturn(0);
414: }

416: /*------------------------------------------------------------*/
420: int TaoCreate_GPCG(TAO_SOLVER tao)
421: {
422:   TAO_GPCG *gpcg;
423:   int      info;

425:   TaoFunctionBegin;

427:   info = TaoNew(TAO_GPCG,&gpcg); CHKERRQ(info);
428:   info = PetscLogObjectMemory(tao,sizeof(TAO_GPCG)); CHKERRQ(info);

430:   info=TaoSetTaoSolveRoutine(tao,TaoSolve_GPCG,(void*)gpcg); CHKERRQ(info);
431:   info=TaoSetTaoSetUpDownRoutines(tao,TaoSetUp_GPCG,TaoSetDown_GPCG); CHKERRQ(info);
432:   info=TaoSetTaoOptionsRoutine(tao,TaoSetFromOptions_GPCG); CHKERRQ(info);
433:   info=TaoSetTaoViewRoutine(tao,TaoView_GPCG); CHKERRQ(info);
434:   info=TaoSetTaoDualVariablesRoutine(tao,TaoGetDualVariables_GPCG); CHKERRQ(info);

436:   info = TaoSetMaximumIterates(tao,500); CHKERRQ(info);
437:   info = TaoSetMaximumFunctionEvaluations(tao,100000); CHKERRQ(info);
438:   info = TaoSetTolerances(tao,1e-12,1e-12,0,0); CHKERRQ(info);

440:   /* Initialize pointers and variables */
441:   gpcg->n=0;
442:   gpcg->maxgpits = 8;
443:   gpcg->pg_ftol = 0.1;

445:   gpcg->gp_iterates=0; /* Cumulative number */
446:   gpcg->total_gp_its = 0;
447:  
448:   /* Initialize pointers and variables */
449:   gpcg->n_bind=0;
450:   gpcg->n_free = 0;
451:   gpcg->n_upper=0;
452:   gpcg->n_lower=0;

454:   //  info = TaoCreateProjectedLineSearch(tao); CHKERRQ(info);
455:   info = TaoGPCGCreateLineSearch(tao); CHKERRQ(info);

457:   TaoFunctionReturn(0);
458: }