Actual source code: hyppilut.c

  1: /*

  3: */

 5:  #include src/ksp/pc/pcimpl.h
  7: #include "HYPRE.h"
  8: #include "IJ_mv.h"
  9: #include "HYPRE_parcsr_ls.h"

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

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

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

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

 35:   /* options for pilut */
 36:   int                factorrowsize;

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

 48:   /* options for euclid */
 49:   PetscTruth         bjilu;
 50:   int                levels;

 52:   /* options for euclid and BoomerAMG */
 53:   PetscTruth         printstatistics;

 55:   /* options for BoomerAMG */
 56:   int                maxlevels;
 57:   double             strongthreshold;
 58:   double             maxrowsum;
 59:   int                *gridsweeps;
 60:   int                coarsentype;
 61:   int                measuretype;
 62:   int                *relaxtype;
 63:   int                **gridrelaxpoints;
 64: } PC_HYPRE;


 69: static PetscErrorCode PCSetUp_HYPRE(PC pc)
 70: {
 71:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
 72:   PetscErrorCode     ierr;
 73:   HYPRE_ParCSRMatrix hmat;
 74:   HYPRE_ParVector    bv,xv;
 75:   int                hierr;

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

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

110: static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
111: {
112:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
113:   PetscErrorCode     ierr;
114:   HYPRE_ParCSRMatrix hmat;
115:   PetscScalar        *bv,*xv;
116:   HYPRE_ParVector    jbv,jxv;
117:   PetscScalar        *sbv,*sxv;
118:   PetscScalar        zero=0.0;
119:   int                hierr;

122:   if (!jac->applyrichardson) {VecSet(&zero,x);}
123:   VecGetArray(b,&bv);
124:   VecGetArray(x,&xv);
125:   HYPREReplacePointer(jac->b,bv,sbv);
126:   HYPREReplacePointer(jac->x,xv,sxv);

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

144: static PetscErrorCode PCDestroy_HYPRE(PC pc)
145: {
146:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;

150:   HYPRE_IJMatrixDestroy(jac->ij);
151:   HYPRE_IJVectorDestroy(jac->b);
152:   HYPRE_IJVectorDestroy(jac->x);
153:   MPI_Comm_free(&(jac->comm_hypre));
154:   (*jac->destroy)(jac->hsolver);
155:   PetscFree(jac);
156:   return(0);
157: }

159: /* --------------------------------------------------------------------------------------------*/
162: static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc)
163: {
164:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
166:   PetscTruth     flag;

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

188: static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
189: {
190:   PC_HYPRE    *jac = (PC_HYPRE*)pc->data;
192:   PetscTruth  iascii;

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

217: /* --------------------------------------------------------------------------------------------*/
220: static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc)
221: {
222:   PC_HYPRE  *jac = (PC_HYPRE*)pc->data;
224:   PetscTruth flag;
225:   char       *args[2];

228:   jac->bjilu              = PETSC_FALSE;
229:   jac->levels             = 1;

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

259: static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
260: {
261:   PC_HYPRE    *jac = (PC_HYPRE*)pc->data;
263:   PetscTruth  iascii;

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

277: /* --------------------------------------------------------------------------------------------*/

279: static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout"};
280: static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
281: static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","","Gauss-Seidel/Jacobi","","","symmetric-Gauss-Seidel/Jacobi",
282:                                             "","","Gaussian-elimination"};
285: static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
286: {
287:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
289:   int            n = 4,i,indx;
290:   PetscTruth     flg;

293:   jac->maxlevels       = 25;
294:   jac->maxiter         = 1;
295:   jac->tol             = 1.e-7;
296:   jac->strongthreshold = .25;
297:   jac->maxrowsum       = .9;
298:   jac->coarsentype     = 6;
299:   jac->measuretype     = 0;
300:   jac->applyrichardson = PETSC_FALSE;

302:   /* this is terrible; HYPRE frees this array so we have to malloc it */
303:   jac->gridsweeps    = (int*)malloc(4*sizeof(int));
304:   jac->gridsweeps[0] = jac->gridsweeps[1] = jac->gridsweeps[2] = 2;
305:   jac->gridsweeps[3] = 1;

307:   jac->relaxtype     = (int*)malloc(4*sizeof(int));
308:   jac->relaxtype[0]  = jac->relaxtype[1] = jac->relaxtype[2] = 3;
309:   jac->relaxtype[3]  = 9;

311:   PetscOptionsHead("HYPRE BoomerAMG Options");
312:     PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);
313:     if (flg) {
314:       if (jac->maxlevels < 2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
315:     }
316:     HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);
317:     PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used","None",jac->maxiter,&jac->maxiter,&flg);
318:     if (flg) {
319:       if (jac->maxiter < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
320:     }
321:     HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);
322:     PetscOptionsScalar("-pc_hypre_boomeramg_tol","Convergence tolerance","None",jac->tol,&jac->tol,&flg);
323:     if (flg) {
324:       if (jac->tol < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %g must be great than or equal zero",jac->tol);
325:     }
326:     HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);
327:     PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);
328:     if (flg) {
329:       if (jac->strongthreshold < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %g must be great than or equal zero",jac->strongthreshold);
330:     }
331:     HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);
332:     PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);
333:     if (flg) {
334:       if (jac->maxrowsum < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %g must be greater than zero",jac->maxrowsum);
335:       if (jac->maxrowsum > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %g must be less than or equal one",jac->maxrowsum);
336:     }
337:     HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);
338: 
339:     n = 4;
340:     PetscOptionsIntArray("-pc_hypre_boomeramg_grid_sweeps","Grid sweeps for fine,down,up,coarse","None",jac->gridsweeps,&n,&flg);
341:     if (flg) {
342:       if (n == 1) {
343:         jac->gridsweeps[1] = jac->gridsweeps[2] = jac->gridsweeps[0];
344:         jac->gridsweeps[3] = 1;
345:         n = 4;
346:       }
347:       if (n != 4) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"You must provide either 1 or 4 values seperated by commas, you provided %d",n);
348:     }
349:     HYPRE_BoomerAMGSetNumGridSweeps(jac->hsolver,jac->gridsweeps);

351:     /*
352:          Suggested by QUANDALLE Philippe <Philippe.QUANDALLE@ifp.fr>

354:         gridrelaxpoints[i][j] are for i=0,1,2,3 (fine,down,up,coarse) and j=sweep number
355:         0 indicates smooth all points
356:         1 indicates smooth coarse points
357:        -1 indicates smooth fine points

359:          Here when j=1 it first smooths all the coarse points, then all the fine points.
360:     */
361:     jac->gridrelaxpoints    = (int**)malloc(4*sizeof(int*));
362:     if(jac->gridsweeps[0]>0) jac->gridrelaxpoints[0] = (int*)malloc(jac->gridsweeps[0]*sizeof(int));
363:     if(jac->gridsweeps[1]>0) jac->gridrelaxpoints[1] = (int*)malloc(jac->gridsweeps[1]*sizeof(int));
364:     if(jac->gridsweeps[2]>0) jac->gridrelaxpoints[2] = (int*)malloc(jac->gridsweeps[2]*sizeof(int));
365:     if(jac->gridsweeps[3]>0) jac->gridrelaxpoints[3] = (int*)malloc(jac->gridsweeps[3]*sizeof(int));

367:     PetscOptionsLogical("-pc_hypre_boomeramg_sweep_all","Sweep all points","None",PETSC_FALSE,&flg,0);
368:     if(jac->gridsweeps[0] == 1) jac->gridrelaxpoints[0][0] = 0;
369:     else if(jac->gridsweeps[0] == 2) {
370:       if (flg) {
371:         jac->gridrelaxpoints[0][0] = 0; jac->gridrelaxpoints[0][1] = 0;
372:       } else {
373:         jac->gridrelaxpoints[0][0] = 1; jac->gridrelaxpoints[0][1] = -1;
374:       }
375:     } else if (jac->gridsweeps[0] > 2) {
376:       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Grid sweeps can only be 0, 1, or 2");
377:     }
378: 
379:     if(jac->gridsweeps[1] == 1) jac->gridrelaxpoints[1][0] = 0;
380:     else if(jac->gridsweeps[1] == 2) {
381:       if (flg) {
382:         jac->gridrelaxpoints[1][0] = 0;  jac->gridrelaxpoints[1][1] = 0;
383:       } else {
384:         jac->gridrelaxpoints[1][0] = 1; jac->gridrelaxpoints[1][1] = -1;
385:       }
386:     } else if (jac->gridsweeps[1] > 2) {
387:       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Grid sweeps can only be 0, 1, or 2");
388:     }

390:     if(jac->gridsweeps[2] == 1) jac->gridrelaxpoints[2][0] = 0;
391:     else if(jac->gridsweeps[2] == 2) {
392:       if (flg) {
393:         jac->gridrelaxpoints[2][0] = 0; jac->gridrelaxpoints[2][1] = 0;
394:       } else {
395:         jac->gridrelaxpoints[2][0] = -1; jac->gridrelaxpoints[2][1] = 1;
396:       }
397:     } else if (jac->gridsweeps[2] > 2) {
398:       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Grid sweeps can only be 0, 1, or 2");
399:     }

401:     for (i=0; i<jac->gridsweeps[3]; i++) {
402:       jac->gridrelaxpoints[3][i] = 0;
403:     }
404:     //HYPRE_BoomerAMGSetGridRelaxPoints(jac->hsolver,jac->gridrelaxpoints);

406:     PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);
407:     if (flg) {
408:       jac->measuretype = indx;
409:     }
410:     HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);
411:     PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,7,HYPREBoomerAMGCoarsenType[6],&indx,&flg);
412:     if (flg) {
413:       jac->coarsentype = indx;
414:     }
415:     HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);
416:     PetscOptionsEList("-pc_hypre_boomeramg_relax_type","Relax type","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);
417:     if (flg) {
418:       jac->relaxtype[0] = jac->relaxtype[1] = jac->relaxtype[2] = indx;
419:     }
420:     PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);
421:     if (flg) {
422:       jac->relaxtype[3] = indx;
423:     }
424:     HYPRE_BoomerAMGSetGridRelaxType(jac->hsolver,jac->relaxtype);
425:     PetscOptionsLogical("-pc_hypre_boomeramg_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);
426:     if (jac->printstatistics) {
427:       HYPRE_BoomerAMGSetPrintLevel(jac->hsolver,3);
428:       HYPRE_BoomerAMGSetDebugFlag(jac->hsolver,3);
429:     }
430:   PetscOptionsTail();
431:   return(0);
432: }

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

442:   PetscLogInfo(pc,"PCApplyRichardson_hypre_BoomerAMG: Warning, convergence critera ignored, using %D iterations\n",its);
443:   HYPRE_BoomerAMGSetMaxIter(jac->hsolver,its);
444:   HYPRE_BoomerAMGSetTol(jac->hsolver,rtol);
445:   jac->applyrichardson = PETSC_TRUE;
446:   PCApply_HYPRE(pc,b,y);
447:   jac->applyrichardson = PETSC_FALSE;
448:   HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);
449:   HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);
450:   return(0);
451: }


456: static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
457: {
458:   PC_HYPRE    *jac = (PC_HYPRE*)pc->data;
460:   PetscTruth  iascii;

463:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
464:   if (iascii) {
465:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");
466:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);
467:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations %d\n",jac->maxiter);
468:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance %g\n",jac->tol);
469:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %g\n",jac->strongthreshold);
470:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %g\n",jac->maxrowsum);
471:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on fine grid %d\n",jac->gridsweeps[0]);
472:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps down         %d\n",jac->gridsweeps[1]);
473:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps up           %d\n",jac->gridsweeps[2]);
474:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on coarse    %d\n",jac->gridsweeps[3]);

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

481:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type    %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);
482:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type    %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);
483:   }
484:   return(0);
485: }

487: /* --------------------------------------------------------------------------------------------*/
490: static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc)
491: {
492:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
494:   int            indx;
495:   PetscTruth     flag;
496:   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};

499:   jac->nlevels     = 1;
500:   jac->threshhold  = .1;
501:   jac->filter      = .1;
502:   jac->loadbal     = 0;
503:   if (PetscLogPrintInfo) {
504:     jac->logging     = (int) PETSC_TRUE;
505:   } else {
506:     jac->logging     = (int) PETSC_FALSE;
507:   }
508:   jac->ruse = (int) PETSC_FALSE;
509:   jac->symt   = 0;

511:   PetscOptionsHead("HYPRE ParaSails Options");
512:     PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);
513:     PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,0);
514:     HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);

516:     PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,0);
517:     HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);

519:     PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,0);
520:     HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);

522:     PetscOptionsLogical("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,0);
523:     HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);

525:     PetscOptionsLogical("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,0);
526:     HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);

528:     PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);
529:     if (flag) {
530:       jac->symt = indx;
531:     }
532:     HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);

534:   PetscOptionsTail();
535:   return(0);
536: }

540: static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
541: {
542:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
544:   PetscTruth     iascii;
545:   const char     *symt;

548:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
549:   if (iascii) {
550:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");
551:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);
552:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %g\n",jac->threshhold);
553:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %g\n",jac->filter);
554:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %g\n",jac->loadbal);
555:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",jac->ruse ? "true" : "false");
556:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",jac->logging ? "true" : "false");
557:     if (!jac->symt) {
558:       symt = "nonsymmetric matrix and preconditioner";
559:     } else if (jac->symt == 1) {
560:       symt = "SPD matrix and preconditioner";
561:     } else if (jac->symt == 2) {
562:       symt = "nonsymmetric matrix but SPD preconditioner";
563:     } else {
564:       SETERRQ1(PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
565:     }
566:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);
567:   }
568:   return(0);
569: }
570: /* ---------------------------------------------------------------------------------*/

575: PetscErrorCode PCHYPRESetType_HYPRE(PC pc,const char name[])
576: {
577:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
579:   PetscTruth     flag;

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

586:   pc->ops->setup          = PCSetUp_HYPRE;
587:   pc->ops->apply          = PCApply_HYPRE;
588:   pc->ops->destroy        = PCDestroy_HYPRE;

590:   jac->maxiter            = PETSC_DEFAULT;
591:   jac->tol                = PETSC_DEFAULT;
592:   jac->printstatistics    = PetscLogPrintInfo;

594:   PetscStrcmp("pilut",name,&flag);
595:   if (flag) {
596:     HYPRE_ParCSRPilutCreate(jac->comm_hypre,&jac->hsolver);
597:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
598:     pc->ops->view           = PCView_HYPRE_Pilut;
599:     jac->destroy            = HYPRE_ParCSRPilutDestroy;
600:     jac->setup              = HYPRE_ParCSRPilutSetup;
601:     jac->solve              = HYPRE_ParCSRPilutSolve;
602:     jac->factorrowsize      = PETSC_DEFAULT;
603:     return(0);
604:   }
605:   PetscStrcmp("parasails",name,&flag);
606:   if (flag) {
607:     HYPRE_ParaSailsCreate(jac->comm_hypre,&jac->hsolver);
608:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
609:     pc->ops->view           = PCView_HYPRE_ParaSails;
610:     jac->destroy            = HYPRE_ParaSailsDestroy;
611:     jac->setup              = HYPRE_ParaSailsSetup;
612:     jac->solve              = HYPRE_ParaSailsSolve;
613:     return(0);
614:   }
615:   PetscStrcmp("euclid",name,&flag);
616:   if (flag) {
617:     HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver);
618:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
619:     pc->ops->view           = PCView_HYPRE_Euclid;
620:     jac->destroy            = HYPRE_EuclidDestroy;
621:     jac->setup              = HYPRE_EuclidSetup;
622:     jac->solve              = HYPRE_EuclidSolve;
623:     return(0);
624:   }
625:   PetscStrcmp("boomeramg",name,&flag);
626:   if (flag) {
627:     HYPRE_BoomerAMGCreate(&jac->hsolver);
628:     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
629:     pc->ops->view            = PCView_HYPRE_BoomerAMG;
630:     jac->destroy             = HYPRE_BoomerAMGDestroy;
631:     jac->setup               = HYPRE_BoomerAMGSetup;
632:     jac->solve               = HYPRE_BoomerAMGSolve;
633:     pc->ops->applyrichardson = PCApplyRichardson_BoomerAMG;
634:     return(0);
635:   }
636:   SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
637:   return(0);
638: }

641: /*
642:     It only gets here if the HYPRE type has not been set before the call to 
643:    ...SetFromOptions() which actually is most of the time
644: */
647: static PetscErrorCode PCSetFromOptions_HYPRE(PC pc)
648: {
650:   int            indx;
651:   const char     *type[] = {"pilut","parasails","boomeramg","euclid"};
652:   PetscTruth     flg;

655:   PetscOptionsHead("HYPRE preconditioner options");
656:     PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"pilut",&indx,&flg);
657:     if (PetscOptionsPublishCount) {   /* force the default if it was not yet set and user did not set with option */
658:       if (!flg && !pc->ops->apply) {
659:         flg   = PETSC_TRUE;
660:         indx = 0;
661:       }
662:     }
663:     if (flg) {
664:       PCHYPRESetType_HYPRE(pc,type[indx]);
665:     }
666:     if (pc->ops->setfromoptions) {
667:       pc->ops->setfromoptions(pc);
668:     }
669:   PetscOptionsTail();
670:   return(0);
671: }

675: /*@C
676:      PCHYPRESetType - Sets which hypre preconditioner you wish to use

678:    Input Parameters:
679: +     pc - the preconditioner context
680: -     name - either  pilut, parasails, boomeramg, euclid

682:    Options Database Keys:
683:    -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
684:  
685:    Level: intermediate

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

690: @*/
691: PetscErrorCode PCHYPRESetType(PC pc,const char name[])
692: {
693:   PetscErrorCode ierr,(*f)(PC,const char[]);

698:   PetscObjectQueryFunction((PetscObject)pc,"PCHYPRESetType_C",(void (**)(void))&f);
699:   if (f) {
700:     (*f)(pc,name);
701:   }
702:   return(0);
703: }

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

708:    Options Database Keys:
709: +   -pc_hypre_type - One of pilut, parasails, boomerAMG, euclid
710: -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
711:           preconditioner
712:  
713:    Level: intermediate

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

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

722: M*/

727: PetscErrorCode PCCreate_HYPRE(PC pc)
728: {
729:   PC_HYPRE       *jac;

733:   PetscNew(PC_HYPRE,&jac);
734:   PetscMemzero(jac,sizeof(PC_HYPRE));
735:   pc->data                 = jac;
736:   pc->ops->setfromoptions  = PCSetFromOptions_HYPRE;
737:   /* Com_dup for hypre */
738:   MPI_Comm_dup(pc->comm,&(jac->comm_hypre));
739:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);
740:   return(0);
741: }