Actual source code: hyppilut.c

  1: #define PETSCKSP_DLL

  3: /*

  5: */

 7:  #include private/pcimpl.h
  9: #include "HYPRE.h"
 10: #include "IJ_mv.h"
 11: #include "parcsr_ls.h"

 14: EXTERN PetscErrorCode MatHYPRE_IJMatrixCreate(Mat,HYPRE_IJMatrix*);
 15: EXTERN PetscErrorCode MatHYPRE_IJMatrixCopy(Mat,HYPRE_IJMatrix);
 16: EXTERN PetscErrorCode VecHYPRE_IJVectorCreate(Vec,HYPRE_IJVector*);

 18: /* 
 19:    Private context (data structure) for the  preconditioner.  
 20: */
 21: typedef struct {
 22:   HYPRE_Solver       hsolver;
 23:   HYPRE_IJMatrix     ij;
 24:   HYPRE_IJVector     b,x;

 26:   PetscErrorCode     (*destroy)(HYPRE_Solver);
 27:   PetscErrorCode     (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
 28:   PetscErrorCode     (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
 29: 
 30:   MPI_Comm           comm_hypre;

 32:   /* options for pilut and BoomerAMG*/
 33:   int                maxiter;
 34:   double             tol;
 35:   PetscTruth         applyrichardson;

 37:   /* options for pilut */
 38:   int                factorrowsize;

 40:   /* options for parasails */
 41:   int                nlevels;
 42:   double             threshhold;
 43:   double             filter;
 44:   int                sym;
 45:   double             loadbal;
 46:   int                logging;
 47:   int                ruse;
 48:   int                symt;

 50:   /* options for euclid */
 51:   PetscTruth         bjilu;
 52:   int                levels;

 54:   /* options for euclid and BoomerAMG */
 55:   PetscTruth         printstatistics;

 57:   /* options for BoomerAMG */
 58:   int                maxlevels;
 59:   double             strongthreshold;
 60:   double             maxrowsum;
 61:   int                gridsweeps[4];
 62:   int                coarsentype;
 63:   int                measuretype;
 64:   int                relaxtype[4];
 65:   double             relaxweight;
 66:   double             outerrelaxweight;
 67:   int                relaxorder;
 68:   int                **gridrelaxpoints;
 69:   double             truncfactor;
 70: } PC_HYPRE;


 75: static PetscErrorCode PCSetUp_HYPRE(PC pc)
 76: {
 77:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
 78:   PetscErrorCode     ierr;
 79:   HYPRE_ParCSRMatrix hmat;
 80:   HYPRE_ParVector    bv,xv;
 81:   int                hierr;

 84:   if (!jac->ij) { /* create the matrix the first time through */
 85:     MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);
 86:   }
 87:   if (!jac->b) {
 88:     Vec vec;
 89:     MatGetVecs(pc->pmat,&vec,0);
 90:     VecHYPRE_IJVectorCreate(vec,&jac->b);
 91:     VecHYPRE_IJVectorCreate(vec,&jac->x);
 92:     VecDestroy(vec);
 93:   }
 94:   MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);
 95:   HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);
 96:   HYPRE_IJVectorGetObject(jac->b,(void**)&bv);
 97:   HYPRE_IJVectorGetObject(jac->x,(void**)&xv);
 98:   h(*jac->setup)(jac->hsolver,hmat,bv,xv);
 99:   if (hierr) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE setup, error code %d",hierr);
100:   return(0);
101: }

103: /*
104:     Replaces the address where the HYPRE vector points to its data with the address of
105:   PETSc's data. Saves the old address so it can be reset when we are finished with it.
106:   Allows use to get the data into a HYPRE vector without the cost of memcopies 
107: */
108: #define HYPREReplacePointer(b,newvalue,savedvalue) {\
109:    hypre_ParVector *par_vector   = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));\
110:    hypre_Vector    *local_vector = hypre_ParVectorLocalVector(par_vector);\
111:    savedvalue         = local_vector->data;\
112:    local_vector->data = newvalue;}

116: static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
117: {
118:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
119:   PetscErrorCode     ierr;
120:   HYPRE_ParCSRMatrix hmat;
121:   PetscScalar        *bv,*xv;
122:   HYPRE_ParVector    jbv,jxv;
123:   PetscScalar        *sbv,*sxv;
124:   int                hierr;

127:   if (!jac->applyrichardson) {VecSet(x,0.0);}
128:   VecGetArray(b,&bv);
129:   VecGetArray(x,&xv);
130:   HYPREReplacePointer(jac->b,bv,sbv);
131:   HYPREReplacePointer(jac->x,xv,sxv);

133:   HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);
134:   HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);
135:   HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);
136:   h(*jac->solve)(jac->hsolver,hmat,jbv,jxv);
137:   /* error code of 1 in boomerAMG merely means convergence not achieved */
138:   if (hierr && (hierr != 1 || jac->solve != HYPRE_BoomerAMGSolve)) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
139:   if (hierr) hypre__global_error = 0;
140: 
141:   HYPREReplacePointer(jac->b,sbv,bv);
142:   HYPREReplacePointer(jac->x,sxv,xv);
143:   VecRestoreArray(x,&xv);
144:   VecRestoreArray(b,&bv);
145:   return(0);
146: }

150: static PetscErrorCode PCDestroy_HYPRE(PC pc)
151: {
152:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;

156:   if (jac->ij) { HYPRE_IJMatrixDestroy(jac->ij); }
157:   if (jac->b) { HYPRE_IJVectorDestroy(jac->b); }
158:   if (jac->x) { HYPRE_IJVectorDestroy(jac->x); }
159:   (*jac->destroy)(jac->hsolver);
160:   MPI_Comm_free(&(jac->comm_hypre));
161:   PetscFree(jac);
162:   return(0);
163: }

165: /* --------------------------------------------------------------------------------------------*/
168: static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc)
169: {
170:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
172:   PetscTruth     flag;

175:   PetscOptionsHead("HYPRE Pilut Options");
176:   PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);
177:   if (flag) {
178:     HYPRE_ParCSRPilutSetMaxIter(jac->hsolver,jac->maxiter);
179:   }
180:   PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);
181:   if (flag) {
182:     HYPRE_ParCSRPilutSetDropTolerance(jac->hsolver,jac->tol);
183:   }
184:   PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);
185:   if (flag) {
186:     HYPRE_ParCSRPilutSetFactorRowSize(jac->hsolver,jac->factorrowsize);
187:   }
188:   PetscOptionsTail();
189:   return(0);
190: }

194: static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
195: {
196:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
198:   PetscTruth     iascii;

201:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
202:   if (iascii) {
203:     PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");
204:     if (jac->maxiter != PETSC_DEFAULT) {
205:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);
206:     } else {
207:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");
208:     }
209:     if (jac->tol != PETSC_DEFAULT) {
210:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %G\n",jac->tol);
211:     } else {
212:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");
213:     }
214:     if (jac->factorrowsize != PETSC_DEFAULT) {
215:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);
216:     } else {
217:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");
218:     }
219:   }
220:   return(0);
221: }

223: /* --------------------------------------------------------------------------------------------*/
226: static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc)
227: {
228:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
230:   PetscTruth     flag;
231:   char           *args[2];

234:   PetscOptionsHead("HYPRE Euclid Options");
235:   PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);
236:   if (flag) {
237:     char levels[16];
238:     if (jac->levels < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels);
239:     sprintf(levels,"%d",jac->levels);
240:     args[0] = (char*)"-level"; args[1] = levels;
241:     HYPRE_EuclidSetParams(jac->hsolver,2,args);
242:   }
243:   PetscOptionsTruth("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);
244:   if (jac->bjilu) {
245:     args[0] =(char*) "-bj"; args[1] = (char*)"1";
246:     HYPRE_EuclidSetParams(jac->hsolver,2,args);
247:   }
248: 
249:   PetscOptionsTruth("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);
250:   if (jac->printstatistics) {
251:     args[0] = (char*)"-eu_stats"; args[1] = (char*)"1";
252:     HYPRE_EuclidSetParams(jac->hsolver,2,args);
253:     args[0] = (char*)"-eu_mem"; args[1] = (char*)"1";
254:     HYPRE_EuclidSetParams(jac->hsolver,2,args);
255:   }
256:   PetscOptionsTail();
257:   return(0);
258: }

262: static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
263: {
264:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
266:   PetscTruth     iascii;

269:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
270:   if (iascii) {
271:     PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n");
272:     PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: number of levels %d\n",jac->levels);
273:     if (jac->bjilu) {
274:       PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");
275:     }
276:   }
277:   return(0);
278: }

280: /* --------------------------------------------------------------------------------------------*/

282: static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout"};
283: static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
284: static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","","SOR/Jacobi","backward-SOR/Jacobi","","symmetric-SOR/Jacobi",
285:                                                   "","","Gaussian-elimination"};
288: static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
289: {
290:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
292:   int            n,indx;
293:   PetscTruth     flg, tmp_truth;
294:   double         tmpdbl, twodbl[2];

297:   PetscOptionsHead("HYPRE BoomerAMG Options");
298:   PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);
299:   if (flg) {
300:     if (jac->maxlevels < 2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
301:     HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);
302:   }
303:   PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);
304:   if (flg) {
305:     if (jac->maxiter < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
306:     HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);
307:   }
308:   PetscOptionsScalar("-pc_hypre_boomeramg_tol","Convergence tolerance PER hypre call","None",jac->tol,&jac->tol,&flg);
309:   if (flg) {
310:     if (jac->tol < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be great than or equal zero",jac->tol);
311:     HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);
312:   }

314:   PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor","None",jac->truncfactor,&jac->truncfactor,&flg);
315:   if (flg) {
316:     if (jac->truncfactor < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor);
317:     HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);
318:   }

320:   PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);
321:   if (flg) {
322:     if (jac->strongthreshold < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold);
323:     HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);
324:   }
325:   PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);
326:   if (flg) {
327:     if (jac->maxrowsum < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum);
328:     if (jac->maxrowsum > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum);
329:     HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);
330:   }

332:   /* Grid sweeps */
333:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_all","Number of sweeps for all grid levels (fine, up, and down)","None", jac->gridsweeps[0], &indx ,&flg);
334:   if (flg) {
335:     HYPRE_BoomerAMGSetNumSweeps(jac->hsolver,indx);
336:     /* modify the jac structure so we can view the updated options with PC_View */
337:     jac->gridsweeps[0] = indx;
338:     jac->gridsweeps[1] = jac->gridsweeps[2] = jac->gridsweeps[0];
339:     jac->gridsweeps[3] = 1;  /*The coarse level is not affected by this function - hypre code sets to 1*/
340:   }
341:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_fine","Number of sweeps for the fine level","None", jac->gridsweeps[0], &indx ,&flg);
342:   if (flg) {
343:     HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 0);
344:     jac->gridsweeps[0] = indx;
345:   }
346:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[2], &indx ,&flg);
347:   if (flg) {
348:     HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 1);
349:     jac->gridsweeps[1] = indx;
350:   }
351:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None", jac->gridsweeps[1], &indx ,&flg);
352:   if (flg) {
353:     HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 2);
354:     jac->gridsweeps[2] = indx;
355:   }
356:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None", jac->gridsweeps[3], &indx ,&flg);
357:   if (flg) {
358:     HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 3);
359:     jac->gridsweeps[3] = indx;
360:   }

362:   /* Relax type */
363:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_all","Relax type for fine, up, and down cycles (coarse level set to gaussian elimination)","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);
364:   if (flg) {
365:     jac->relaxtype[0] = jac->relaxtype[1] = jac->relaxtype[2] = indx;
366:     jac->relaxtype[3] = 9; /* hypre code sets coarse grid to 9 (G.E.)*/
367:     hypre_BoomerAMGSetRelaxType(jac->hsolver, indx);
368:   }
369:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_fine","Relax type on fine grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);
370:   if (flg) {
371:     jac->relaxtype[0] = indx;
372:     HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 0);
373:   }
374:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);
375:   if (flg) {
376:     jac->relaxtype[1] = indx;
377:     HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 1);
378:   }
379:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);
380:   if (flg) {
381:     jac->relaxtype[2] = indx;
382:     HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 2);
383:   }
384:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);
385:   if (flg) {
386:     jac->relaxtype[3] = indx;
387:     HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 3);
388:   }

390:   /* Relaxation Weight */
391:   PetscOptionsReal("-pc_hypre_boomeramg_relax_weight_all","Relaxation weight for all levels (0 = hypre estimates, -k = determined with k CG steps)","None",jac->relaxweight, &tmpdbl ,&flg);
392:   if (flg) {
393:     hypre_BoomerAMGSetRelaxWt( jac->hsolver, tmpdbl);
394:     jac->relaxweight = tmpdbl;
395:   }

397:   n=2;
398:   twodbl[0] = twodbl[1] = 1.0;
399:   PetscOptionsRealArray("-pc_hypre_boomeramg_relax_weight_level","Set the relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);
400:   if (flg) {
401:     if (n == 2) {
402:       indx =  (int)PetscAbsReal(twodbl[1]);
403:       hypre_BoomerAMGSetLevelRelaxWt( jac->hsolver, twodbl[0], indx);
404:     } else {
405:       SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Relax weight level: you must provide 2 values separated by a comma (and no space), you provided %d",n);
406:     }
407:   }

409:   /* Outer relaxation Weight */
410:   PetscOptionsReal("-pc_hypre_boomeramg_outer_relax_weight_all","Outer relaxation weight for all levels ( -k = determined with k CG steps)","None",jac->outerrelaxweight, &tmpdbl ,&flg);
411:   if (flg) {
412:     hypre_BoomerAMGSetOuterWt( jac->hsolver, tmpdbl);
413:     jac->outerrelaxweight = tmpdbl;
414:   }

416:   n=2;
417:   twodbl[0] = twodbl[1] = 1.0;
418:   PetscOptionsRealArray("-pc_hypre_boomeramg_outer_relax_weight_level","Set the outer relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);
419:   if (flg) {
420:     if (n == 2) {
421:       indx =  (int)PetscAbsReal(twodbl[1]);
422:       hypre_BoomerAMGSetLevelOuterWt( jac->hsolver, twodbl[0], indx);
423:     } else {
424:       SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Relax weight outer level: You must provide 2 values separated by a comma (and no space), you provided %d",n);
425:     }
426:   }

428:   /* the Relax Order */
429:   /* PetscOptionsName("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", &flg); */
430:   PetscOptionsTruth( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);

432:   if (flg) {
433:     jac->relaxorder = 0;
434:     hypre_BoomerAMGSetRelaxOrder( jac->hsolver, jac->relaxorder);
435:   }
436:   PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);
437:   if (flg) {
438:     jac->measuretype = indx;
439:     HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);
440:   }
441:   PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,7,HYPREBoomerAMGCoarsenType[6],&indx,&flg);
442:   if (flg) {
443:     jac->coarsentype = indx;
444:     HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);
445:   }
446:   PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);
447:   if (flg) {
448:     int level=3;
449:     jac->printstatistics = PETSC_TRUE;
450:     PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);
451:     HYPRE_BoomerAMGSetPrintLevel(jac->hsolver,level);
452:     HYPRE_BoomerAMGSetDebugFlag(jac->hsolver,level);
453:   }
454:   PetscOptionsTail();
455:   return(0);
456: }

460: static PetscErrorCode PCApplyRichardson_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
461: {
462:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;

466:   HYPRE_BoomerAMGSetMaxIter(jac->hsolver,its*jac->maxiter);
467:   HYPRE_BoomerAMGSetTol(jac->hsolver,PetscMin(rtol,jac->tol));
468:   jac->applyrichardson = PETSC_TRUE;
469:   PCApply_HYPRE(pc,b,y);
470:   jac->applyrichardson = PETSC_FALSE;
471:   HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);
472:   HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);
473:   return(0);
474: }


479: static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
480: {
481:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
483:   PetscTruth     iascii;

486:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
487:   if (iascii) {
488:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");
489:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);
490:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations PER hypre call%d\n",jac->maxiter);
491:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance PER hypre call%G\n",jac->tol);
492:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);
493:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);

495:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on fine grid %d\n",jac->gridsweeps[0]);
496:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps down         %d\n",jac->gridsweeps[1]);
497:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps up           %d\n",jac->gridsweeps[2]);
498:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on coarse    %d\n",jac->gridsweeps[3]);

500:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on fine grid  %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);
501:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);
502:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);
503:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[3]]);

505:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax weight  (all)      %G\n",jac->relaxweight);
506:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);

508:     if (jac->relaxorder) {
509:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using CF-relaxation\n");
510:     } else {
511:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using CF-relaxation\n");
512:     }
513:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);
514:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);
515:   }
516:   return(0);
517: }

519: /* --------------------------------------------------------------------------------------------*/
522: static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc)
523: {
524:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
526:   int            indx;
527:   PetscTruth     flag;
528:   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};

531:   PetscOptionsHead("HYPRE ParaSails Options");
532:   PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);
533:   PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);
534:   if (flag) {
535:     HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);
536:   }

538:   PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);
539:   if (flag) {
540:     HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);
541:   }

543:   PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);
544:   if (flag) {
545:     HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);
546:   }

548:   PetscOptionsTruth("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,&flag);
549:   if (flag) {
550:     HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);
551:   }

553:   PetscOptionsTruth("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,&flag);
554:   if (flag) {
555:     HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);
556:   }

558:   PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);
559:   if (flag) {
560:     jac->symt = indx;
561:     HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);
562:   }

564:   PetscOptionsTail();
565:   return(0);
566: }

570: static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
571: {
572:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
574:   PetscTruth     iascii;
575:   const char     *symt = 0;;

578:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
579:   if (iascii) {
580:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");
581:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);
582:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %G\n",jac->threshhold);
583:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %G\n",jac->filter);
584:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %G\n",jac->loadbal);
585:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscTruths[jac->ruse]);
586:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscTruths[jac->logging]);
587:     if (!jac->symt) {
588:       symt = "nonsymmetric matrix and preconditioner";
589:     } else if (jac->symt == 1) {
590:       symt = "SPD matrix and preconditioner";
591:     } else if (jac->symt == 2) {
592:       symt = "nonsymmetric matrix but SPD preconditioner";
593:     } else {
594:       SETERRQ1(PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
595:     }
596:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);
597:   }
598:   return(0);
599: }
600: /* ---------------------------------------------------------------------------------*/

605: PetscErrorCode  PCHYPRESetType_HYPRE(PC pc,const char name[])
606: {
607:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
609:   PetscTruth     flag;
610:   PetscInt       bs;

613:   if (pc->ops->setup) {
614:     SETERRQ(PETSC_ERR_ORDER,"Cannot set the HYPRE preconditioner type once it has been set");
615:   }

617:   pc->ops->setup          = PCSetUp_HYPRE;
618:   pc->ops->apply          = PCApply_HYPRE;
619:   pc->ops->destroy        = PCDestroy_HYPRE;

621:   jac->maxiter            = PETSC_DEFAULT;
622:   jac->tol                = PETSC_DEFAULT;
623:   jac->printstatistics    = PetscLogPrintInfo;

625:   PetscStrcmp("pilut",name,&flag);
626:   if (flag) {
627:     HYPRE_ParCSRPilutCreate(jac->comm_hypre,&jac->hsolver);
628:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
629:     pc->ops->view           = PCView_HYPRE_Pilut;
630:     jac->destroy            = HYPRE_ParCSRPilutDestroy;
631:     jac->setup              = HYPRE_ParCSRPilutSetup;
632:     jac->solve              = HYPRE_ParCSRPilutSolve;
633:     jac->factorrowsize      = PETSC_DEFAULT;
634:     return(0);
635:   }
636:   PetscStrcmp("parasails",name,&flag);
637:   if (flag) {
638:     HYPRE_ParaSailsCreate(jac->comm_hypre,&jac->hsolver);
639:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
640:     pc->ops->view           = PCView_HYPRE_ParaSails;
641:     jac->destroy            = HYPRE_ParaSailsDestroy;
642:     jac->setup              = HYPRE_ParaSailsSetup;
643:     jac->solve              = HYPRE_ParaSailsSolve;
644:     /* initialize */
645:     jac->nlevels     = 1;
646:     jac->threshhold  = .1;
647:     jac->filter      = .1;
648:     jac->loadbal     = 0;
649:     if (PetscLogPrintInfo) {
650:       jac->logging     = (int) PETSC_TRUE;
651:     } else {
652:       jac->logging     = (int) PETSC_FALSE;
653:     }
654:     jac->ruse = (int) PETSC_FALSE;
655:     jac->symt   = 0;
656:     HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);
657:     HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);
658:     HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);
659:     HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);
660:     HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);
661:     HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);
662:     return(0);
663:   }
664:   PetscStrcmp("euclid",name,&flag);
665:   if (flag) {
666:     HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver);
667:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
668:     pc->ops->view           = PCView_HYPRE_Euclid;
669:     jac->destroy            = HYPRE_EuclidDestroy;
670:     jac->setup              = HYPRE_EuclidSetup;
671:     jac->solve              = HYPRE_EuclidSolve;
672:     /* initialization */
673:     jac->bjilu              = PETSC_FALSE;
674:     jac->levels             = 1;
675:     return(0);
676:   }
677:   PetscStrcmp("boomeramg",name,&flag);
678:   if (flag) {
679:     HYPRE_BoomerAMGCreate(&jac->hsolver);
680:     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
681:     pc->ops->view            = PCView_HYPRE_BoomerAMG;
682:     jac->destroy             = HYPRE_BoomerAMGDestroy;
683:     jac->setup               = HYPRE_BoomerAMGSetup;
684:     jac->solve               = HYPRE_BoomerAMGSolve;
685:     pc->ops->applyrichardson = PCApplyRichardson_BoomerAMG;
686:     /* these defaults match the hypre defaults */
687:     jac->maxlevels        = 25;
688:     jac->maxiter          = 1;
689:     jac->tol              = 1.e-7;
690:     jac->strongthreshold  = .25;
691:     jac->maxrowsum        = .9;
692:     jac->coarsentype      = 6;
693:     jac->measuretype      = 0;
694:     jac->applyrichardson  = PETSC_FALSE;
695:     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = jac->gridsweeps[3]  = 1;
696:     jac->relaxtype[0]     = jac->relaxtype[1]  = jac->relaxtype[2]  = 3;
697:     jac->relaxtype[3]     = 9;
698:     jac->relaxweight      = 1.0;
699:     jac->outerrelaxweight = 1.0;
700:     jac->relaxorder       = 1;
701:     HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);
702:     HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);
703:     HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);
704:     HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);
705:     HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);
706:     HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);
707:     HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);
708:     HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);
709:     MatGetBlockSize(pc->pmat,&bs);
710:     HYPRE_BoomerAMGSetNumFunctions(jac->hsolver,bs);
711:     return(0);
712:   }
713:   SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
714:   return(0);
715: }

718: /*
719:     It only gets here if the HYPRE type has not been set before the call to 
720:    ...SetFromOptions() which actually is most of the time
721: */
724: static PetscErrorCode PCSetFromOptions_HYPRE(PC pc)
725: {
727:   int            indx;
728:   const char     *type[] = {"pilut","parasails","boomeramg","euclid"};
729:   PetscTruth     flg;

732:   PetscOptionsHead("HYPRE preconditioner options");
733:   PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"pilut",&indx,&flg);
734:   if (PetscOptionsPublishCount) {   /* force the default if it was not yet set and user did not set with option */
735:     if (!flg && !pc->ops->apply) {
736:       flg   = PETSC_TRUE;
737:       indx = 0;
738:     }
739:   }
740:   if (flg) {
741:     PCHYPRESetType_HYPRE(pc,type[indx]);
742:   }
743:   if (pc->ops->setfromoptions) {
744:     pc->ops->setfromoptions(pc);
745:   }
746:   PetscOptionsTail();
747:   return(0);
748: }

752: /*@C
753:      PCHYPRESetType - Sets which hypre preconditioner you wish to use

755:    Input Parameters:
756: +     pc - the preconditioner context
757: -     name - either  pilut, parasails, boomeramg, euclid

759:    Options Database Keys:
760:    -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
761:  
762:    Level: intermediate

764: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
765:            PCHYPRE

767: @*/
768: PetscErrorCode  PCHYPRESetType(PC pc,const char name[])
769: {
770:   PetscErrorCode ierr,(*f)(PC,const char[]);

775:   PetscObjectQueryFunction((PetscObject)pc,"PCHYPRESetType_C",(void (**)(void))&f);
776:   if (f) {
777:     (*f)(pc,name);
778:   }
779:   return(0);
780: }

782: /*MC
783:      PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre

785:    Options Database Keys:
786: +   -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
787: -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
788:           preconditioner
789:  
790:    Level: intermediate

792:    Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()),
793:           the many hypre options can ONLY be set via the options database (e.g. the command line
794:           or with PetscOptionsSetValue(), there are no functions to set them)

796:           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations
797:           (V-cycles) that boomeramg does EACH time it is called. So for example, if it is set to 2 then 
798:           2-V-cycles are being used to define the preconditioner. -ksp_max_iter and -ksp_rtol STILL determine
799:           the total number of iterations and tolerance for the Krylov solver. For example, if 
800:           -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10 then AT MOST twenty V-cycles of boomeramg
801:           will be called.


804:           If you wish to use boomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly
805:           and use -ksp_max_it to control the number of V-cycles.
806:           (see the PETSc FAQ.html at the PETSc website under the Documentation tab).

808: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
809:            PCHYPRESetType()

811: M*/

816: PetscErrorCode  PCCreate_HYPRE(PC pc)
817: {
818:   PC_HYPRE       *jac;

822:   PetscNew(PC_HYPRE,&jac);
823:   pc->data                 = jac;
824:   pc->ops->setfromoptions  = PCSetFromOptions_HYPRE;
825:   /* Com_dup for hypre */
826:   MPI_Comm_dup(pc->comm,&(jac->comm_hypre));
827:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);
828:   return(0);
829: }