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