Actual source code: daapp_mgrid.c

  1: #include "taodaapplication.h"     /* taodaapplication.h */
  2: //#include "taoapp.h"
  3: #include "src/petsctao/include/taopetsc.h"

  5: #include "../../application/taoapp/taoapp_petsc.h"
  6: #include "../interface/daapp_impl.h"
  7: #include "multigridmat.h"

  9: #include "petscmg.h"

 11: class TaoMultiGridApplication: public TaoPetscApplication {

 13:  protected:

 15:  public:

 17:   TaoMultiGridApplication(MPI_Comm);
 18:   TaoMultiGridApplication();

 20:   int coarselevel;
 21:   int EvaluateHessian(TaoVec *xx, TaoMat *HH);

 23: };


 28: TaoMultiGridApplication :: TaoMultiGridApplication(MPI_Comm mpicomm) : TaoPetscApplication(mpicomm) {
 29:   Mat M=0;
 30:   this->coarselevel=0;
 31:   this->taoh=new TaoMultiGridMatPetsc(M);
 32:   this->taoj=new TaoMultiGridMatPetsc(M);
 33:   return;
 34: }

 38: TaoMultiGridApplication::TaoMultiGridApplication() : TaoPetscApplication(PETSC_COMM_WORLD){
 39:   
 40:   Mat M=0;
 41:   this->coarselevel=0;
 42:   this->taoh=new TaoMultiGridMatPetsc(M);
 43:   this->taoj=new TaoMultiGridMatPetsc(M);
 44:   return;
 45: }


 50: int TaoMultiGridApplication::EvaluateHessian(TaoVec *xx, TaoMat *HH)
 51: {
 52:   TaoVecPetsc *px = dynamic_cast <TaoVecPetsc *> (xx);
 53:   TaoMatPetsc *pH = dynamic_cast <TaoMatPetsc *> (HH);

 55:   int     info,i,clevel,cclevel;
 56:   Vec X = px->GetVec();
 57:   Mat H, HPre;
 58:   MatStructure pflag;
 59:   DA_APPLICATION daapp;

 62:   info=TaoAppDAApp(this->papp,&daapp); CHKERRQ(info);

 64:   PetscStackPush("TAO User DA Hessian");
 65:   info=pH->GetMatrix(&H,&HPre,&pflag);CHKERRQ(info);
 66:   info=DAAppGetCurrentLevel(this->papp,&clevel); CHKERRQ(info);
 67:   cclevel=PetscMin(clevel,this->coarselevel);
 68:   cclevel=PetscMax(cclevel,0);

 70:   for (i=clevel; i>=cclevel; i--){
 71:     daapp->currentlevel=i;
 72:     info=TaoAppSetColoring(this->papp,daapp->grid[i].coloring); CHKERRQ(info);
 73:     info = PetscInfo1(daapp,"TAO hessian evaluation at DA_APPLICATION object, level %d.\n",i); CHKERRQ(info);
 74:     info = TaoAppComputeHessian(this->papp, X, &H, &HPre, &pflag); CHKERRQ(info);
 75:     if (i==clevel){  
 76:       info = pH->SetMatrix(H,HPre,pflag);CHKERRQ(info);
 77:     }

 79:     if (i>cclevel){
 80:       info=MatMultTranspose(daapp->grid[i].Interpolate,X,daapp->grid[i-1].R);CHKERRQ(info);
 81:       info=VecPointwiseMult(daapp->grid[i-1].R,daapp->grid[i].CScale,daapp->grid[i-1].R);CHKERRQ(info);
 82:       X=daapp->grid[i-1].R;
 83:       H=daapp->grid[i-1].H;
 84:       HPre=H;
 85:     }
 86:   }
 87:   daapp->currentlevel=clevel;
 88:   info=TaoAppSetColoring(this->papp,daapp->grid[clevel].coloring); CHKERRQ(info);
 89:   PetscStackPop;

 91:   return(0);
 92: }

 94: int DAAppSetMultigridKSP(GridCtx *, int, KSP);
 95: int DAAppSetupMultigridMonitor(TAO_APPLICATION, DA, int, void *);
 96: int TaoAppGetMultiGridApplication(TAO_APPLICATION, TaoMultiGridApplication *);
 97: static char TaoPetscPCMGDAAppXXX[] = "TaoPetscPCMGDAApp";
 98: static char TaoPetscAppXXX[] = "TaoPetscApp";


103: static int DAAppDestroyTaoAppXXX(void*ctx){
104:   TaoMultiGridApplication *mctx=(TaoMultiGridApplication*)ctx;
106:   if (mctx){ delete mctx;}
107:   return(0);
108: }

112: int TaoAppGetMultiGridApplication(TAO_APPLICATION daapplication, TaoMultiGridApplication **mgdaapp){
113:   int ii,clevel,info;
114:   Mat H,HH;
115:   MatStructure         pflag=SAME_NONZERO_PATTERN;
116:   DA_APPLICATION       daapp;
117:   TaoMultiGridApplication*   daapp2;
118:   TaoPetscApplication  *pppappp;

122:   info=TaoAppQueryForObject(daapplication,TaoPetscPCMGDAAppXXX,(void**)&daapp2); CHKERRQ(info);
123:   if (daapp2==0){
124:     daapp2=new TaoMultiGridApplication(PETSC_COMM_WORLD);
125:     daapp2->papp=daapplication;
126:     info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);

128:     info=TaoAppQueryForObject(daapplication,TaoPetscAppXXX,(void**)&pppappp); CHKERRQ(info);
129:     if (pppappp && pppappp->taoh && 0==1){
130:       info=pppappp->taoh->GetMatrix(&H,&HH,&pflag); CHKERRQ(info);
131:       info=daapp2->taoh->SetMatrix(H,HH,pflag);CHKERRQ(info);
132:     } else {
133:       info=DAAppGetCurrentLevel(daapplication,&clevel); CHKERRQ(info);
134:       info=daapp2->taoh->SetMatrix(daapp->grid[clevel].H,daapp->grid[clevel].H,pflag);CHKERRQ(info);
135:     }
136:     if (pppappp){
137:       daapp2->tao=pppappp->tao;
138:     }
139:     info=TaoAppAddObject(daapplication,TaoPetscPCMGDAAppXXX,(void*)daapp2,&ii); CHKERRQ(info);
140:     info=TaoAppQueryRemoveObject(daapplication,TaoPetscAppXXX); CHKERRQ(info);
141:     info=TaoAppAddObject(daapplication,TaoPetscAppXXX,(void*)daapp2,&ii); CHKERRQ(info);
142:     info=TaoAppSetDestroyRoutine(daapplication,DAAppDestroyTaoAppXXX, (void*)daapp2); CHKERRQ(info);
143:   }
144:   *mgdaapp=daapp2;
145:   return(0);
146: }

148: static int UsePCMGGG=0;
151: /*@
152:   DAAppUseMultigrid - Set the preconditioner for the linear solver to be an algebraic multigrid.

154:   Collective on TAO_APPLICATION
155:   
156:   Input Parameters:
157: +  daapplication - the DA Application object
158: -  coarselevel - the coarsest grid to be used in the multigrid preconditioner. (Grid 0 is the coarsest grid.

160:    Level: intermediate

162:    Options Database Key:
163: +  -tao_da_multigrid - use multigrid linear solver
164: -  -ksp_view - view the linear solver

166: Note:
167:   This function should be called after DAAppSetHessianRoutine();

169: Note:
170:   This function should be called before each optimization solver as part of the DAAppMonitor

172: Note:
173:   Multigrid functionality is still under developement for good performance.

175: .seealso: TaoAppGetKSP(), DAAppSetupMultigrid()

177:    Options Database Key:
178: .  -tao_da_multigrid

180: .keywords: Linear Solver, Multigrid, DA, KSP

182: @*/
183: int DAAppUseMultigrid(TAO_APPLICATION daapplication, int coarselevel){
184:   int info;
186:   UsePCMGGG=coarselevel;
187:   info=DAAppSetBeforeMonitor(daapplication,DAAppSetupMultigridMonitor,0); 
188:   CHKERRQ(info);
189:   return(0);
190: }


195: int DAAppSetupMultigridMonitor(TAO_APPLICATION daapplication, DA da, int level, void *ctx){
196:   int info;
198:   info = DAAppSetupMultigrid(daapplication,UsePCMGGG); CHKERRQ(info);
199:   return(0);
200: }


205: /*@
206:    DAAppSetupMultigrid - Sets up the multigrid preconditioner.

208:    Collective on TAO_APPLICATION

210:    Input Parameters:
211: +  daapplication - the TAO_APPLICATION context
212: -  coarselevel - the coarsest grid that should be used in the multigrid (>=0)

214:    Note:
215:    Usually the coarselevel is set to 0;

217:    Level: intermediate

219: .seealso:  DAAppUseMultigrid(), TaoAppSetMonitor()
220: @*/
221: int DAAppSetupMultigrid(TAO_APPLICATION daapplication, int coarselevel){
222:   int nn,level,info;
223:   TaoMultiGridMatPetsc *MMM;
224:   DA_APPLICATION daapp;
225:   TaoMultiGridApplication *mgdaapp;
226:   TAO_SOLVER tao;
227:   KSP  ksp;
228:  
230:   info=DAAppGetCurrentLevel(daapplication,&level); CHKERRQ(info);
231:   info=TaoAppDAApp(daapplication,&daapp); CHKERRQ(info);
232:   if (daapp->grid[level].mgrid==0 &&  coarselevel<level){
233:     info=TaoAppGetMultiGridApplication(daapplication,&mgdaapp); CHKERRQ(info);
234:     mgdaapp->coarselevel=PetscMax(0,coarselevel);
235:     nn=level-coarselevel+1;
236:     if ( mgdaapp->taoh ){
237:       MMM=(TaoMultiGridMatPetsc*)(mgdaapp->taoh);
238:       info=MMM->SetUpMultiGrid(daapp->grid+coarselevel,nn); CHKERRQ(info);
239:     } else if ( mgdaapp->taoj ){
240:       MMM=(TaoMultiGridMatPetsc*)(mgdaapp->taoj);
241:       info=MMM->SetUpMultiGrid(daapp->grid+coarselevel,nn); CHKERRQ(info);
242:       mgdaapp->taoj=MMM;
243:     }

245:     info=TaoAppGetTaoSolver(daapplication,&tao); CHKERRQ(info);
246:     info=TaoSetupApplicationSolver(daapplication,tao); CHKERRQ(info);
247:     info=TaoAppGetKSP(daapplication,&ksp); CHKERRQ(info);
248:     if (ksp){
250:       info = DAAppSetMultigridKSP(daapp->grid+coarselevel,nn,ksp);CHKERRQ(info);
251:     }
252:     daapp->grid[level].mgrid=1;
253:   }
254:   return(0);
255: }


258:   
261: int DAAppSetMultigridKSP(GridCtx* dagrid, int nlevels, KSP ksp){
262:   int i,info;
263:   PC pc,pc2;
264:   KSP sksp;
265:   KSPType ksptype;
266:   PCType pctype;
267:   PCMGType form;

270:   if (ksp==0 || nlevels<=1){
271:     return(0);
272:   }

275:   info = PetscInfo(ksp,"Set multigrid precondition in optimization solver.\n"); CHKERRQ(info);
276:   
277:   info = KSPGetType(ksp,&ksptype); CHKERRQ(info);
278:   
279:   info = KSPGetPC(ksp,&pc); CHKERRQ(info);
280:   info = PCGetType(pc,&pctype); CHKERRQ(info);
281:   info = PCSetType(pc,PCMG); CHKERRQ(info);
282:   
283:   info = PCMGSetLevels(pc,nlevels,PETSC_NULL); CHKERRQ(info);

285:   for (i=0; i<nlevels; i++) {
286:     
287:     info = PCMGGetSmoother(pc,i,&sksp); CHKERRQ(info);
288:     info = KSPSetType(sksp,ksptype); CHKERRQ(info);

290:     info = KSPGetPC(sksp,&pc2); CHKERRQ(info);
291:     info = PCSetType(pc2,PCJACOBI); CHKERRQ(info);

293:     form=PC_MG_MULTIPLICATIVE;
294:     info=PCMGSetType(pc,form);CHKERRQ(info);

296:     info = KSPSetOperators(sksp,dagrid[i].H,dagrid[i].H,SAME_NONZERO_PATTERN); CHKERRQ(info);
297:     
298:     if (i<nlevels-1){
299:       info = PCMGSetX(pc,i,dagrid[i].X); CHKERRQ(info);
300:       info = PCMGSetRhs(pc,i,dagrid[i].RHS); CHKERRQ(info);      
301:     }
302:     if (i>0){
303:       info = PCMGSetR(pc,i,dagrid[i].R); CHKERRQ(info);
304:     }
305:     info = PCMGSetResidual(pc,i,PCMGDefaultResidual,dagrid[i].H); CHKERRQ(info);
306:   }
307:   for (i=1; i<nlevels; i++) {
308:     info = PCMGSetInterpolation(pc,i,dagrid[i].Interpolate); CHKERRQ(info);
309:     info = PCMGSetRestriction(pc,i,dagrid[i].Interpolate); CHKERRQ(info);
310:   }
311:   
312:   return(0);
313: }

317: /* @
318:   DAAppSetQuadraticObjective - identify the objective function as quadratic or not.

320:   Collective on TAO_APPLICATION

322:   Input Parameters:
323: +  daapplication - the DA Application object
324: -  flag - indicates whether the objective is quadratic or not

326:    Level: intermediate

328:    Options Database Key:
329: .  -tao_da_quadratic - use multigrid linear solver

331:    Note:
332:    The default value is PETSC_FALSE.

334:    Note:
335:    If quadratic, consider setting the flag used in KSP to SAME_PRECONDITIONER

337:    Note:
338:    If quadratic, this routine reduces the number of Hessian evaluations done when using
339:    the multigrid preconditioner.

341: .seealso: DAAppSetMatStructure(),  DAAppUseMultigrid()


344: .keywords: DA, Objective Function

346: @ */
347: int DAAppSetQuadraticObjective(TAO_APPLICATION daapplication, PetscTruth pflag){
348:   int info;
349:   TaoMultiGridApplication *mgdaapp;

352:   info=TaoAppGetMultiGridApplication(daapplication,&mgdaapp); CHKERRQ(info);
353:   //  mgdaapp->IsQuadratic=pflag;
354:   if (pflag==PETSC_TRUE){
355:     info = PetscInfo(daapplication,"User labels the objective function as quadratic.\n"); CHKERRQ(info);
356:     info = PetscInfo(daapplication,"Set KSP MatStructure to SAME_PRECONDITIONER.\n"); CHKERRQ(info);
357:     info=DAAppSetMatStructure(daapplication, SAME_PRECONDITIONER);CHKERRQ(info);
358:   } else {
359:     info = PetscInfo(daapplication,"User labels the objective function as NOT quadratic.\n"); CHKERRQ(info);
360:     info = PetscInfo(daapplication,"Set KSP MatStructure to SAME_NONZERO_PATTERN.\n"); CHKERRQ(info);
361:     info=DAAppSetMatStructure(daapplication, SAME_NONZERO_PATTERN);CHKERRQ(info);
362:   }
363:   return(0);
364: }


369: /* @
370:   DAAppSetMultiGridOptions - Sets various multigrid options to be used in this application
371:   and the TAO solver.

373:    Collective on TAO_APPLICATION

375:    Input Parameters:
376: .  daapplication - the DA Application object

378:    Level: beginner

380: .keywords:  options

382: .seealso: TaoDAAppSetOptions();

384: @ */
385: int DAAppSetMultiGridOptions(TAO_APPLICATION daapplication){
386:   int info,nlevel;

388:   PetscTruth flg1=PETSC_FALSE;

391:   info = PetscInfo(daapplication,"TaoDAAppSetMultiGridOptions(): Reading command line for options\n"); CHKERRQ(info);

393:   flg1=PETSC_FALSE;
394:   info = PetscOptionsInt("-tao_da_multigrid","use multigrid linear solver","DAAppUseMultigrid",
395:                          PETSC_FALSE,&nlevel,&flg1);CHKERRQ(info);
396:   if (flg1) {
397:     info=DAAppUseMultigrid(daapplication,nlevel);CHKERRQ(info);
398:     info = PetscInfo(daapplication,"TaoDAAppSetOptions: Use Multigrid linear solver \n"); CHKERRQ(info);
399:   }
400:   /*
401:   flg1=PETSC_FALSE,flg2=PETSC_FALSE;
402:   info = PetscOptionsTruth("-tao_da_quadratic","Identify the objective function as quadratic",
403:                              "DAAppSetQuadraticObjective",PETSC_FALSE,&flg2,&flg1);CHKERRQ(info);
404:   if (flg1) {
405:     info = DAAppSetQuadraticObjective(daapplication, flg2);CHKERRQ(info); 
406:   }
407:   */

409:   return(0);
410: }