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