xref: /petsc/src/ksp/pc/impls/hypre/hypre.c (revision 71f87433c2d5ca8a2055e7bac4e3ff8159d57028)
1 #define PETSCKSP_DLL
2 
3 /*
4    Provides an interface to the LLNL package hypre
5 */
6 #include "private/pcimpl.h"          /*I "petscpc.h" I*/
7 EXTERN_C_BEGIN
8 #include "HYPRE.h"
9 #include "IJ_mv.h"
10 #include "parcsr_ls.h"
11 EXTERN_C_END
12 
13 EXTERN PetscErrorCode MatHYPRE_IJMatrixCreate(Mat,HYPRE_IJMatrix*);
14 EXTERN PetscErrorCode MatHYPRE_IJMatrixCopy(Mat,HYPRE_IJMatrix);
15 EXTERN PetscErrorCode VecHYPRE_IJVectorCreate(Vec,HYPRE_IJVector*);
16 
17 /*
18    Private context (data structure) for the  preconditioner.
19 */
20 typedef struct {
21   HYPRE_Solver       hsolver;
22   HYPRE_IJMatrix     ij;
23   HYPRE_IJVector     b,x;
24 
25   PetscErrorCode     (*destroy)(HYPRE_Solver);
26   PetscErrorCode     (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
27   PetscErrorCode     (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
28 
29   MPI_Comm           comm_hypre;
30   char              *hypre_type;
31 
32   /* options for Pilut and BoomerAMG*/
33   int                maxiter;
34   double             tol;
35 
36   /* options for Pilut */
37   int                factorrowsize;
38 
39   /* options for ParaSails */
40   int                nlevels;
41   double             threshhold;
42   double             filter;
43   int                sym;
44   double             loadbal;
45   int                logging;
46   int                ruse;
47   int                symt;
48 
49   /* options for Euclid */
50   PetscTruth         bjilu;
51   int                levels;
52 
53   /* options for Euclid and BoomerAMG */
54   PetscTruth         printstatistics;
55 
56   /* options for BoomerAMG */
57   int                cycletype;
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   double             truncfactor;
69   PetscTruth         applyrichardson;
70 
71 } PC_HYPRE;
72 
73 
74 #undef __FUNCT__
75 #define __FUNCT__ "PCSetUp_HYPRE"
76 static PetscErrorCode PCSetUp_HYPRE(PC pc)
77 {
78   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
79   PetscErrorCode     ierr;
80   HYPRE_ParCSRMatrix hmat;
81   HYPRE_ParVector    bv,xv;
82   PetscInt           bs;
83   int                hierr;
84 
85   PetscFunctionBegin;
86   if (!jac->hypre_type) {
87     ierr = PCHYPRESetType(pc,"pilut");CHKERRQ(ierr);
88   }
89 #if 1
90   if (!pc->setupcalled) {
91     /* create the matrix and vectors the first time through */
92     Vec x,b;
93     ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr);
94     ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr);
95     ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr);
96     ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr);
97     ierr = VecDestroy(x);CHKERRQ(ierr);
98     ierr = VecDestroy(b);CHKERRQ(ierr);
99   } else if (pc->flag != SAME_NONZERO_PATTERN) {
100     /* rebuild the matrix from scratch */
101     ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr);
102     ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr);
103   }
104 #else
105   if (!jac->ij) { /* create the matrix the first time through */
106     ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr);
107   }
108   if (!jac->b) { /* create the vectors the first time through */
109     Vec x,b;
110     ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr);
111     ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr);
112     ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr);
113     ierr = VecDestroy(x);CHKERRQ(ierr);
114     ierr = VecDestroy(b);CHKERRQ(ierr);
115   }
116 #endif
117   /* special case for BoomerAMG */
118   if (jac->setup == HYPRE_BoomerAMGSetup) {
119     ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr);
120     if (bs > 1) {
121       ierr = HYPRE_BoomerAMGSetNumFunctions(jac->hsolver,bs);CHKERRQ(ierr);
122     }
123   };
124   ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr);
125   ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr);
126   ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&bv);CHKERRQ(ierr);
127   ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&xv);CHKERRQ(ierr);
128   hierr = (*jac->setup)(jac->hsolver,hmat,bv,xv);
129   if (hierr) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE setup, error code %d",hierr);
130   PetscFunctionReturn(0);
131 }
132 
133 /*
134     Replaces the address where the HYPRE vector points to its data with the address of
135   PETSc's data. Saves the old address so it can be reset when we are finished with it.
136   Allows use to get the data into a HYPRE vector without the cost of memcopies
137 */
138 #define HYPREReplacePointer(b,newvalue,savedvalue) {\
139    hypre_ParVector *par_vector   = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));\
140    hypre_Vector    *local_vector = hypre_ParVectorLocalVector(par_vector);\
141    savedvalue         = local_vector->data;\
142    local_vector->data = newvalue;}
143 
144 #undef __FUNCT__
145 #define __FUNCT__ "PCApply_HYPRE"
146 static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
147 {
148   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
149   PetscErrorCode     ierr;
150   HYPRE_ParCSRMatrix hmat;
151   PetscScalar        *bv,*xv;
152   HYPRE_ParVector    jbv,jxv;
153   PetscScalar        *sbv,*sxv;
154   int                hierr;
155 
156   PetscFunctionBegin;
157   if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);}
158   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
159   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
160   HYPREReplacePointer(jac->b,bv,sbv);
161   HYPREReplacePointer(jac->x,xv,sxv);
162 
163   ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr);
164   ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr);
165   ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr);
166   hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv);
167   /* error code of 1 in BoomerAMG merely means convergence not achieved */
168   if (hierr && (hierr != 1 || jac->solve != HYPRE_BoomerAMGSolve)) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
169   if (hierr) hypre__global_error = 0;
170 
171   HYPREReplacePointer(jac->b,sbv,bv);
172   HYPREReplacePointer(jac->x,sxv,xv);
173   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
174   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
175   PetscFunctionReturn(0);
176 }
177 
178 #undef __FUNCT__
179 #define __FUNCT__ "PCDestroy_HYPRE"
180 static PetscErrorCode PCDestroy_HYPRE(PC pc)
181 {
182   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
183   PetscErrorCode ierr;
184 
185   PetscFunctionBegin;
186   if (jac->ij) { ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); }
187   if (jac->b)  { ierr = HYPRE_IJVectorDestroy(jac->b);CHKERRQ(ierr);  }
188   if (jac->x)  { ierr = HYPRE_IJVectorDestroy(jac->x);CHKERRQ(ierr);  }
189   if (jac->destroy) { ierr = (*jac->destroy)(jac->hsolver);CHKERRQ(ierr); }
190   ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr);
191   if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);}
192   ierr = PetscFree(jac);CHKERRQ(ierr);
193 
194   ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr);
195   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);CHKERRQ(ierr);
196   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);CHKERRQ(ierr);
197   PetscFunctionReturn(0);
198 }
199 
200 /* --------------------------------------------------------------------------------------------*/
201 #undef __FUNCT__
202 #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut"
203 static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc)
204 {
205   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
206   PetscErrorCode ierr;
207   PetscTruth     flag;
208 
209   PetscFunctionBegin;
210   ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr);
211   ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr);
212   if (flag) {
213     ierr = HYPRE_ParCSRPilutSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
214   }
215   ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr);
216   if (flag) {
217     ierr = HYPRE_ParCSRPilutSetDropTolerance(jac->hsolver,jac->tol);CHKERRQ(ierr);
218   }
219   ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr);
220   if (flag) {
221     ierr = HYPRE_ParCSRPilutSetFactorRowSize(jac->hsolver,jac->factorrowsize);CHKERRQ(ierr);
222   }
223   ierr = PetscOptionsTail();CHKERRQ(ierr);
224   PetscFunctionReturn(0);
225 }
226 
227 #undef __FUNCT__
228 #define __FUNCT__ "PCView_HYPRE_Pilut"
229 static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
230 {
231   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
232   PetscErrorCode ierr;
233   PetscTruth     iascii;
234 
235   PetscFunctionBegin;
236   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
237   if (iascii) {
238     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");CHKERRQ(ierr);
239     if (jac->maxiter != PETSC_DEFAULT) {
240       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr);
241     } else {
242       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr);
243     }
244     if (jac->tol != PETSC_DEFAULT) {
245       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr);
246     } else {
247       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr);
248     }
249     if (jac->factorrowsize != PETSC_DEFAULT) {
250       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr);
251     } else {
252       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");CHKERRQ(ierr);
253     }
254   }
255   PetscFunctionReturn(0);
256 }
257 
258 /* --------------------------------------------------------------------------------------------*/
259 #undef __FUNCT__
260 #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid"
261 static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc)
262 {
263   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
264   PetscErrorCode ierr;
265   PetscTruth     flag;
266   char           *args[8],levels[16];
267   PetscInt       cnt = 0;
268 
269   PetscFunctionBegin;
270   ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr);
271   ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr);
272   if (flag) {
273     if (jac->levels < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels);
274     sprintf(levels,"%d",jac->levels);
275     args[cnt++] = (char*)"-level"; args[cnt++] = levels;
276   }
277   ierr = PetscOptionsTruth("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);CHKERRQ(ierr);
278   if (jac->bjilu) {
279     args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1";
280   }
281 
282   ierr = PetscOptionsTruth("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);CHKERRQ(ierr);
283   if (jac->printstatistics) {
284     args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1";
285     args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1";
286   }
287   ierr = PetscOptionsTail();CHKERRQ(ierr);
288   if (cnt) {
289     ierr = HYPRE_EuclidSetParams(jac->hsolver,cnt,args);CHKERRQ(ierr);
290   }
291   PetscFunctionReturn(0);
292 }
293 
294 #undef __FUNCT__
295 #define __FUNCT__ "PCView_HYPRE_Euclid"
296 static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
297 {
298   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
299   PetscErrorCode ierr;
300   PetscTruth     iascii;
301 
302   PetscFunctionBegin;
303   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
304   if (iascii) {
305     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n");CHKERRQ(ierr);
306     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr);
307     if (jac->bjilu) {
308       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr);
309     }
310   }
311   PetscFunctionReturn(0);
312 }
313 
314 /* --------------------------------------------------------------------------------------------*/
315 
316 #undef __FUNCT__
317 #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG"
318 static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
319 {
320   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
321   PetscErrorCode     ierr;
322   HYPRE_ParCSRMatrix hmat;
323   PetscScalar        *bv,*xv;
324   HYPRE_ParVector    jbv,jxv;
325   PetscScalar        *sbv,*sxv;
326   int                hierr;
327 
328   PetscFunctionBegin;
329   ierr = VecSet(x,0.0);CHKERRQ(ierr);
330   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
331   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
332   HYPREReplacePointer(jac->b,bv,sbv);
333   HYPREReplacePointer(jac->x,xv,sxv);
334 
335   ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr);
336   ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr);
337   ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr);
338 
339   hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv);
340   /* error code of 1 in BoomerAMG merely means convergence not achieved */
341   if (hierr && (hierr != 1)) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
342   if (hierr) hypre__global_error = 0;
343 
344   HYPREReplacePointer(jac->b,sbv,bv);
345   HYPREReplacePointer(jac->x,sxv,xv);
346   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
347   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
348   PetscFunctionReturn(0);
349 }
350 
351 static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
352 static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout"};
353 static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
354 static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","","SOR/Jacobi","backward-SOR/Jacobi","","symmetric-SOR/Jacobi","","","Gaussian-elimination"};
355 
356 #undef __FUNCT__
357 #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG"
358 static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
359 {
360   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
361   PetscErrorCode ierr;
362   int            n,indx;
363   PetscTruth     flg, tmp_truth;
364   double         tmpdbl, twodbl[2];
365 
366   PetscFunctionBegin;
367   ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr);
368   ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr);
369   if (flg) {
370     jac->cycletype = indx;
371     ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr);
372   }
373   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr);
374   if (flg) {
375     if (jac->maxlevels < 2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
376     ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr);
377   }
378   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr);
379   if (flg) {
380     if (jac->maxiter < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
381     ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
382   }
383   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_tol","Convergence tolerance PER hypre call","None",jac->tol,&jac->tol,&flg);CHKERRQ(ierr);
384   if (flg) {
385     if (jac->tol < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be great than or equal zero",jac->tol);
386     ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
387   }
388 
389   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr);
390   if (flg) {
391     if (jac->truncfactor < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor);
392     ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr);
393   }
394 
395   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr);
396   if (flg) {
397     if (jac->strongthreshold < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold);
398     ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr);
399   }
400   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr);
401   if (flg) {
402     if (jac->maxrowsum < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum);
403     if (jac->maxrowsum > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum);
404     ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr);
405   }
406 
407   /* Grid sweeps */
408   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_all","Number of sweeps for all grid levels (fine, up, and down)","None",jac->gridsweeps[0],&indx,&flg);CHKERRQ(ierr);
409   if (flg) {
410     ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver,indx);CHKERRQ(ierr);
411     /* modify the jac structure so we can view the updated options with PC_View */
412     jac->gridsweeps[0] = indx;
413     jac->gridsweeps[1] = jac->gridsweeps[2] = jac->gridsweeps[0];
414     jac->gridsweeps[3] = 1;  /*The coarse level is not affected by this function - hypre code sets to 1*/
415   }
416   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_fine","Number of sweeps for the fine level","None",jac->gridsweeps[0],&indx,&flg);CHKERRQ(ierr);
417   if (flg) {
418     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 0);CHKERRQ(ierr);
419     jac->gridsweeps[0] = indx;
420   }
421   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr);
422   if (flg) {
423     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 1);CHKERRQ(ierr);
424     jac->gridsweeps[1] = indx;
425   }
426   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr);
427   if (flg) {
428     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 2);CHKERRQ(ierr);
429     jac->gridsweeps[2] = indx;
430   }
431   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[3],&indx,&flg);CHKERRQ(ierr);
432   if (flg) {
433     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 3);CHKERRQ(ierr);
434     jac->gridsweeps[3] = indx;
435   }
436 
437   /* Relax type */
438   ierr = 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);CHKERRQ(ierr);
439   if (flg) {
440     jac->relaxtype[0] = jac->relaxtype[1] = jac->relaxtype[2] = indx;
441     jac->relaxtype[3] = 9; /* hypre code sets coarse grid to 9 (G.E.)*/
442     ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, indx);CHKERRQ(ierr);
443   }
444   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_fine","Relax type on fine grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);CHKERRQ(ierr);
445   if (flg) {
446     jac->relaxtype[0] = indx;
447     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 0);CHKERRQ(ierr);
448   }
449   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);CHKERRQ(ierr);
450   if (flg) {
451     jac->relaxtype[1] = indx;
452     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 1);CHKERRQ(ierr);
453   }
454   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);CHKERRQ(ierr);
455   if (flg) {
456     jac->relaxtype[2] = indx;
457     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 2);CHKERRQ(ierr);
458   }
459   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr);
460   if (flg) {
461     jac->relaxtype[3] = indx;
462     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 3);CHKERRQ(ierr);
463   }
464 
465   /* Relaxation Weight */
466   ierr = 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);CHKERRQ(ierr);
467   if (flg) {
468     ierr = HYPRE_BoomerAMGSetRelaxWt(jac->hsolver,tmpdbl);CHKERRQ(ierr);
469     jac->relaxweight = tmpdbl;
470   }
471 
472   n=2;
473   twodbl[0] = twodbl[1] = 1.0;
474   ierr = PetscOptionsRealArray("-pc_hypre_boomeramg_relax_weight_level","Set the relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);CHKERRQ(ierr);
475   if (flg) {
476     if (n == 2) {
477       indx =  (int)PetscAbsReal(twodbl[1]);
478       ierr = HYPRE_BoomerAMGSetLevelRelaxWt(jac->hsolver,twodbl[0],indx);CHKERRQ(ierr);
479     } else {
480       SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Relax weight level: you must provide 2 values separated by a comma (and no space), you provided %d",n);
481     }
482   }
483 
484   /* Outer relaxation Weight */
485   ierr = 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);CHKERRQ(ierr);
486   if (flg) {
487     ierr = HYPRE_BoomerAMGSetOuterWt( jac->hsolver, tmpdbl);CHKERRQ(ierr);
488     jac->outerrelaxweight = tmpdbl;
489   }
490 
491   n=2;
492   twodbl[0] = twodbl[1] = 1.0;
493   ierr = PetscOptionsRealArray("-pc_hypre_boomeramg_outer_relax_weight_level","Set the outer relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);CHKERRQ(ierr);
494   if (flg) {
495     if (n == 2) {
496       indx =  (int)PetscAbsReal(twodbl[1]);
497       ierr = HYPRE_BoomerAMGSetLevelOuterWt( jac->hsolver, twodbl[0], indx);CHKERRQ(ierr);
498     } else {
499       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);
500     }
501   }
502 
503   /* the Relax Order */
504   /* ierr = PetscOptionsName("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", &flg);CHKERRQ(ierr); */
505   ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
506 
507   if (flg) {
508     jac->relaxorder = 0;
509     ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr);
510   }
511   ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr);
512   if (flg) {
513     jac->measuretype = indx;
514     ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr);
515   }
516   ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,7,HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr);
517   if (flg) {
518     jac->coarsentype = indx;
519     ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr);
520   }
521   ierr = PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);CHKERRQ(ierr);
522   if (flg) {
523     int level=3;
524     jac->printstatistics = PETSC_TRUE;
525     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
526     ierr = HYPRE_BoomerAMGSetPrintLevel(jac->hsolver,level);CHKERRQ(ierr);
527     ierr = HYPRE_BoomerAMGSetDebugFlag(jac->hsolver,level);CHKERRQ(ierr);
528   }
529   ierr = PetscOptionsTail();CHKERRQ(ierr);
530   PetscFunctionReturn(0);
531 }
532 
533 #undef __FUNCT__
534 #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG"
535 static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
536 {
537   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
538   PetscErrorCode ierr;
539 
540   PetscFunctionBegin;
541   ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,its*jac->maxiter);CHKERRQ(ierr);
542   ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,PetscMin(rtol,jac->tol));CHKERRQ(ierr);
543   jac->applyrichardson = PETSC_TRUE;
544   ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr);
545   jac->applyrichardson = PETSC_FALSE;
546   ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
547   ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
548   PetscFunctionReturn(0);
549 }
550 
551 
552 #undef __FUNCT__
553 #define __FUNCT__ "PCView_HYPRE_BoomerAMG"
554 static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
555 {
556   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
557   PetscErrorCode ierr;
558   PetscTruth     iascii;
559 
560   PetscFunctionBegin;
561   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
562   if (iascii) {
563     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr);
564     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr);
565     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr);
566     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr);
567     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr);
568     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr);
569     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr);
570 
571     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on fine grid %d\n",jac->gridsweeps[0]);CHKERRQ(ierr);
572     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps down         %d\n",jac->gridsweeps[1]);CHKERRQ(ierr);
573     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps up           %d\n",jac->gridsweeps[2]);CHKERRQ(ierr);
574     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on coarse    %d\n",jac->gridsweeps[3]);CHKERRQ(ierr);
575 
576     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on fine grid  %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr);
577     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr);
578     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr);
579     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[3]]);CHKERRQ(ierr);
580 
581     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax weight  (all)      %G\n",jac->relaxweight);CHKERRQ(ierr);
582     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr);
583 
584     if (jac->relaxorder) {
585       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr);
586     } else {
587       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr);
588     }
589     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr);
590     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr);
591   }
592   PetscFunctionReturn(0);
593 }
594 
595 /* --------------------------------------------------------------------------------------------*/
596 #undef __FUNCT__
597 #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails"
598 static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc)
599 {
600   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
601   PetscErrorCode ierr;
602   int            indx;
603   PetscTruth     flag;
604   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};
605 
606   PetscFunctionBegin;
607   ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr);
608   ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr);
609   ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr);
610   if (flag) {
611     ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr);
612   }
613 
614   ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr);
615   if (flag) {
616     ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr);
617   }
618 
619   ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr);
620   if (flag) {
621     ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr);
622   }
623 
624   ierr = PetscOptionsTruth("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,&flag);CHKERRQ(ierr);
625   if (flag) {
626     ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr);
627   }
628 
629   ierr = PetscOptionsTruth("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,&flag);CHKERRQ(ierr);
630   if (flag) {
631     ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr);
632   }
633 
634   ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);CHKERRQ(ierr);
635   if (flag) {
636     jac->symt = indx;
637     ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr);
638   }
639 
640   ierr = PetscOptionsTail();CHKERRQ(ierr);
641   PetscFunctionReturn(0);
642 }
643 
644 #undef __FUNCT__
645 #define __FUNCT__ "PCView_HYPRE_ParaSails"
646 static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
647 {
648   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
649   PetscErrorCode ierr;
650   PetscTruth     iascii;
651   const char     *symt = 0;;
652 
653   PetscFunctionBegin;
654   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
655   if (iascii) {
656     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");CHKERRQ(ierr);
657     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr);
658     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr);
659     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr);
660     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr);
661     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscTruths[jac->ruse]);CHKERRQ(ierr);
662     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscTruths[jac->logging]);CHKERRQ(ierr);
663     if (!jac->symt) {
664       symt = "nonsymmetric matrix and preconditioner";
665     } else if (jac->symt == 1) {
666       symt = "SPD matrix and preconditioner";
667     } else if (jac->symt == 2) {
668       symt = "nonsymmetric matrix but SPD preconditioner";
669     } else {
670       SETERRQ1(PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
671     }
672     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr);
673   }
674   PetscFunctionReturn(0);
675 }
676 /* ---------------------------------------------------------------------------------*/
677 
678 EXTERN_C_BEGIN
679 #undef __FUNCT__
680 #define __FUNCT__ "PCHYPREGetType_HYPRE"
681 PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType_HYPRE(PC pc,const char *name[])
682 {
683   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
684 
685   PetscFunctionBegin;
686   *name = jac->hypre_type;
687   PetscFunctionReturn(0);
688 }
689 EXTERN_C_END
690 
691 EXTERN_C_BEGIN
692 #undef __FUNCT__
693 #define __FUNCT__ "PCHYPRESetType_HYPRE"
694 PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType_HYPRE(PC pc,const char name[])
695 {
696   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
697   PetscErrorCode ierr;
698   PetscTruth     flag;
699 
700   PetscFunctionBegin;
701   if (jac->hypre_type) {
702     ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr);
703     if (!flag) SETERRQ(PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
704     PetscFunctionReturn(0);
705   } else {
706     ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr);
707   }
708 
709   jac->maxiter            = PETSC_DEFAULT;
710   jac->tol                = PETSC_DEFAULT;
711   jac->printstatistics    = PetscLogPrintInfo;
712 
713   ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr);
714   if (flag) {
715     ierr                    = HYPRE_ParCSRPilutCreate(jac->comm_hypre,&jac->hsolver);
716     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
717     pc->ops->view           = PCView_HYPRE_Pilut;
718     jac->destroy            = HYPRE_ParCSRPilutDestroy;
719     jac->setup              = HYPRE_ParCSRPilutSetup;
720     jac->solve              = HYPRE_ParCSRPilutSolve;
721     jac->factorrowsize      = PETSC_DEFAULT;
722     PetscFunctionReturn(0);
723   }
724   ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr);
725   if (flag) {
726     ierr                    = HYPRE_ParaSailsCreate(jac->comm_hypre,&jac->hsolver);
727     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
728     pc->ops->view           = PCView_HYPRE_ParaSails;
729     jac->destroy            = HYPRE_ParaSailsDestroy;
730     jac->setup              = HYPRE_ParaSailsSetup;
731     jac->solve              = HYPRE_ParaSailsSolve;
732     /* initialize */
733     jac->nlevels     = 1;
734     jac->threshhold  = .1;
735     jac->filter      = .1;
736     jac->loadbal     = 0;
737     if (PetscLogPrintInfo) {
738       jac->logging     = (int) PETSC_TRUE;
739     } else {
740       jac->logging     = (int) PETSC_FALSE;
741     }
742     jac->ruse = (int) PETSC_FALSE;
743     jac->symt   = 0;
744     ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr);
745     ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr);
746     ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr);
747     ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr);
748     ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr);
749     ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr);
750     PetscFunctionReturn(0);
751   }
752   ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr);
753   if (flag) {
754     ierr                    = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver);
755     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
756     pc->ops->view           = PCView_HYPRE_Euclid;
757     jac->destroy            = HYPRE_EuclidDestroy;
758     jac->setup              = HYPRE_EuclidSetup;
759     jac->solve              = HYPRE_EuclidSolve;
760     /* initialization */
761     jac->bjilu              = PETSC_FALSE;
762     jac->levels             = 1;
763     PetscFunctionReturn(0);
764   }
765   ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr);
766   if (flag) {
767     ierr                     = HYPRE_BoomerAMGCreate(&jac->hsolver);
768     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
769     pc->ops->view            = PCView_HYPRE_BoomerAMG;
770     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
771     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
772     jac->destroy             = HYPRE_BoomerAMGDestroy;
773     jac->setup               = HYPRE_BoomerAMGSetup;
774     jac->solve               = HYPRE_BoomerAMGSolve;
775     jac->applyrichardson     = PETSC_FALSE;
776     /* these defaults match the hypre defaults */
777     jac->cycletype        = 1;
778     jac->maxlevels        = 25;
779     jac->maxiter          = 1;
780     jac->tol              = 1.e-7;
781     jac->truncfactor      = 0.0;
782     jac->strongthreshold  = .25;
783     jac->maxrowsum        = .9;
784     jac->coarsentype      = 6;
785     jac->measuretype      = 0;
786     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = jac->gridsweeps[3]  = 1;
787     jac->relaxtype[0]     = jac->relaxtype[1]  = jac->relaxtype[2]  = 3;
788     jac->relaxtype[3]     = 9;
789     jac->relaxweight      = 1.0;
790     jac->outerrelaxweight = 1.0;
791     jac->relaxorder       = 1;
792     ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr);
793     ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr);
794     ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
795     ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
796     ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr);
797     ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr);
798     ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr);
799     ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr);
800     ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr);
801     ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr);
802     PetscFunctionReturn(0);
803   }
804   ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr);
805   jac->hypre_type = PETSC_NULL;
806   SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
807   PetscFunctionReturn(0);
808 }
809 EXTERN_C_END
810 
811 /*
812     It only gets here if the HYPRE type has not been set before the call to
813    ...SetFromOptions() which actually is most of the time
814 */
815 #undef __FUNCT__
816 #define __FUNCT__ "PCSetFromOptions_HYPRE"
817 static PetscErrorCode PCSetFromOptions_HYPRE(PC pc)
818 {
819   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
820   PetscErrorCode ierr;
821   int            indx;
822   const char     *type[] = {"pilut","parasails","boomeramg","euclid"};
823   PetscTruth     flg;
824 
825   PetscFunctionBegin;
826   ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr);
827   ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"pilut",&indx,&flg);CHKERRQ(ierr);
828   if (PetscOptionsPublishCount) {   /* force the default if it was not yet set and user did not set with option */
829     if (!flg && !jac->hypre_type) {
830       flg   = PETSC_TRUE;
831       indx = 0;
832     }
833   }
834   if (flg) {
835     ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr);
836   }
837   if (pc->ops->setfromoptions) {
838     ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr);
839   }
840   ierr = PetscOptionsTail();CHKERRQ(ierr);
841   PetscFunctionReturn(0);
842 }
843 
844 #undef __FUNCT__
845 #define __FUNCT__ "PCHYPRESetType"
846 /*@C
847      PCHYPRESetType - Sets which hypre preconditioner you wish to use
848 
849    Input Parameters:
850 +     pc - the preconditioner context
851 -     name - either  pilut, parasails, boomeramg, euclid
852 
853    Options Database Keys:
854    -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
855 
856    Level: intermediate
857 
858 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
859            PCHYPRE
860 
861 @*/
862 PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC pc,const char name[])
863 {
864   PetscErrorCode ierr,(*f)(PC,const char[]);
865 
866   PetscFunctionBegin;
867   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
868   PetscValidCharPointer(name,2);
869   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPRESetType_C",(void (**)(void))&f);CHKERRQ(ierr);
870   if (f) {
871     ierr = (*f)(pc,name);CHKERRQ(ierr);
872   }
873   PetscFunctionReturn(0);
874 }
875 
876 #undef __FUNCT__
877 #define __FUNCT__ "PCHYPREGetType"
878 /*@C
879      PCHYPREGetType - Gets which hypre preconditioner you are using
880 
881    Input Parameter:
882 .     pc - the preconditioner context
883 
884    Output Parameter:
885 .     name - either  pilut, parasails, boomeramg, euclid
886 
887    Level: intermediate
888 
889 .seealso:  PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC,
890            PCHYPRE
891 
892 @*/
893 PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC pc,const char *name[])
894 {
895   PetscErrorCode ierr,(*f)(PC,const char*[]);
896 
897   PetscFunctionBegin;
898   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
899   PetscValidPointer(name,2);
900   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPREGetType_C",(void (**)(void))&f);CHKERRQ(ierr);
901   if (f) {
902     ierr = (*f)(pc,name);CHKERRQ(ierr);
903   }
904   PetscFunctionReturn(0);
905 }
906 
907 /*MC
908      PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre
909 
910    Options Database Keys:
911 +   -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
912 -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
913           preconditioner
914 
915    Level: intermediate
916 
917    Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()),
918           the many hypre options can ONLY be set via the options database (e.g. the command line
919           or with PetscOptionsSetValue(), there are no functions to set them)
920 
921           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations
922           (V-cycles) that boomeramg does EACH time it is called. So for example, if it is set to 2 then
923           2-V-cycles are being used to define the preconditioner. -ksp_max_iter and -ksp_rtol STILL determine
924           the total number of iterations and tolerance for the Krylov solver. For example, if
925           -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10 then AT MOST twenty V-cycles of boomeramg
926           will be called.
927 
928           If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly
929           and use -ksp_max_it to control the number of V-cycles.
930           (see the PETSc FAQ.html at the PETSc website under the Documentation tab).
931 
932 	  2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option
933 	  -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L.
934 
935 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
936            PCHYPRESetType()
937 
938 M*/
939 
940 EXTERN_C_BEGIN
941 #undef __FUNCT__
942 #define __FUNCT__ "PCCreate_HYPRE"
943 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_HYPRE(PC pc)
944 {
945   PC_HYPRE       *jac;
946   PetscErrorCode ierr;
947 
948   PetscFunctionBegin;
949   ierr = PetscNew(PC_HYPRE,&jac);CHKERRQ(ierr);
950   ierr = PetscLogObjectMemory(pc,sizeof(PC_HYPRE));CHKERRQ(ierr);
951   pc->data                 = jac;
952   pc->ops->destroy         = PCDestroy_HYPRE;
953   pc->ops->setfromoptions  = PCSetFromOptions_HYPRE;
954   pc->ops->setup           = PCSetUp_HYPRE;
955   pc->ops->apply           = PCApply_HYPRE;
956   jac->comm_hypre          = MPI_COMM_NULL;
957   jac->hypre_type          = PETSC_NULL;
958   /* duplicate communicator for hypre */
959   ierr = MPI_Comm_dup(pc->comm,&(jac->comm_hypre));CHKERRQ(ierr);
960   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr);
961   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr);
962   PetscFunctionReturn(0);
963 }
964 EXTERN_C_END
965