xref: /petsc/src/ksp/pc/impls/hypre/hypre.c (revision dd2cd8d73f22f6c4ff65b606e93d00fe4b08e8e4)
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,"boomeramg");CHKERRQ(ierr);
99   }
100 
101   if (pc->setupcalled) {
102     /* always destroy the old matrix and create a new memory;
103        hope this does not churn the memory too much. The problem
104        is I do not know if it is possible to put the matrix back to
105        its initial state so that we can directly copy the values
106        the second time through. */
107     ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr);
108     jac->ij = 0;
109   }
110 
111   if (!jac->ij) { /* create the matrix the first time through */
112     ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr);
113   }
114   if (!jac->b) { /* create the vectors the first time through */
115     Vec x,b;
116     ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr);
117     ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr);
118     ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr);
119     ierr = VecDestroy(x);CHKERRQ(ierr);
120     ierr = VecDestroy(b);CHKERRQ(ierr);
121   }
122 
123   /* special case for BoomerAMG */
124   if (jac->setup == HYPRE_BoomerAMGSetup) {
125     ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr);
126     if (bs > 1) {
127       ierr = HYPRE_BoomerAMGSetNumFunctions(jac->hsolver,bs);CHKERRQ(ierr);
128     }
129   };
130   ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr);
131   ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr);
132   ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&bv);CHKERRQ(ierr);
133   ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&xv);CHKERRQ(ierr);
134   hierr = (*jac->setup)(jac->hsolver,hmat,bv,xv);
135   if (hierr) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE setup, error code %d",hierr);
136   PetscFunctionReturn(0);
137 }
138 
139 /*
140     Replaces the address where the HYPRE vector points to its data with the address of
141   PETSc's data. Saves the old address so it can be reset when we are finished with it.
142   Allows use to get the data into a HYPRE vector without the cost of memcopies
143 */
144 #define HYPREReplacePointer(b,newvalue,savedvalue) {\
145    hypre_ParVector *par_vector   = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));\
146    hypre_Vector    *local_vector = hypre_ParVectorLocalVector(par_vector);\
147    savedvalue         = local_vector->data;\
148    local_vector->data = newvalue;}
149 
150 #undef __FUNCT__
151 #define __FUNCT__ "PCApply_HYPRE"
152 static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
153 {
154   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
155   PetscErrorCode     ierr;
156   HYPRE_ParCSRMatrix hmat;
157   PetscScalar        *bv,*xv;
158   HYPRE_ParVector    jbv,jxv;
159   PetscScalar        *sbv,*sxv;
160   int                hierr;
161 
162   PetscFunctionBegin;
163   if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);}
164   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
165   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
166   HYPREReplacePointer(jac->b,bv,sbv);
167   HYPREReplacePointer(jac->x,xv,sxv);
168 
169   ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr);
170   ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr);
171   ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr);
172   hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv);
173 
174   /*if (hierr && (hierr != HYPRE_ERROR_CONV || jac->solve != HYPRE_BoomerAMGSolve))SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
175    */
176  /* error code of HYPRE_ERROR_CONV means convergence not achieved - if
177     the tolerance is set to 0.0 (the default), a convergence error will
178     not occur (so we may not want to overide the conv. error here?*/
179  if (hierr && hierr != HYPRE_ERROR_CONV)
180   {
181      SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
182   }
183  if (hierr) hypre__global_error = 0;
184 
185 
186   HYPREReplacePointer(jac->b,sbv,bv);
187   HYPREReplacePointer(jac->x,sxv,xv);
188   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
189   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
190   PetscFunctionReturn(0);
191 }
192 
193 #undef __FUNCT__
194 #define __FUNCT__ "PCDestroy_HYPRE"
195 static PetscErrorCode PCDestroy_HYPRE(PC pc)
196 {
197   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
198   PetscErrorCode ierr;
199 
200   PetscFunctionBegin;
201   if (jac->ij) { ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); }
202   if (jac->b)  { ierr = HYPRE_IJVectorDestroy(jac->b);CHKERRQ(ierr);  }
203   if (jac->x)  { ierr = HYPRE_IJVectorDestroy(jac->x);CHKERRQ(ierr);  }
204   if (jac->destroy) { ierr = (*jac->destroy)(jac->hsolver);CHKERRQ(ierr); }
205   ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr);
206   if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);}
207   ierr = PetscFree(jac);CHKERRQ(ierr);
208 
209   ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr);
210   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);CHKERRQ(ierr);
211   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);CHKERRQ(ierr);
212   PetscFunctionReturn(0);
213 }
214 
215 /* --------------------------------------------------------------------------------------------*/
216 #undef __FUNCT__
217 #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut"
218 static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc)
219 {
220   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
221   PetscErrorCode ierr;
222   PetscTruth     flag;
223 
224   PetscFunctionBegin;
225   ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr);
226   ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr);
227   if (flag) {
228     ierr = HYPRE_ParCSRPilutSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
229   }
230   ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr);
231   if (flag) {
232     ierr = HYPRE_ParCSRPilutSetDropTolerance(jac->hsolver,jac->tol);CHKERRQ(ierr);
233   }
234   ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr);
235   if (flag) {
236     ierr = HYPRE_ParCSRPilutSetFactorRowSize(jac->hsolver,jac->factorrowsize);CHKERRQ(ierr);
237   }
238   ierr = PetscOptionsTail();CHKERRQ(ierr);
239   PetscFunctionReturn(0);
240 }
241 
242 #undef __FUNCT__
243 #define __FUNCT__ "PCView_HYPRE_Pilut"
244 static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
245 {
246   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
247   PetscErrorCode ierr;
248   PetscTruth     iascii;
249 
250   PetscFunctionBegin;
251   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
252   if (iascii) {
253     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");CHKERRQ(ierr);
254     if (jac->maxiter != PETSC_DEFAULT) {
255       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr);
256     } else {
257       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr);
258     }
259     if (jac->tol != PETSC_DEFAULT) {
260       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr);
261     } else {
262       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr);
263     }
264     if (jac->factorrowsize != PETSC_DEFAULT) {
265       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr);
266     } else {
267       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");CHKERRQ(ierr);
268     }
269   }
270   PetscFunctionReturn(0);
271 }
272 
273 /* --------------------------------------------------------------------------------------------*/
274 #undef __FUNCT__
275 #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid"
276 static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc)
277 {
278   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
279   PetscErrorCode ierr;
280   PetscTruth     flag;
281   char           *args[8],levels[16];
282   PetscInt       cnt = 0;
283 
284   PetscFunctionBegin;
285   ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr);
286   ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr);
287   if (flag) {
288     if (jac->levels < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels);
289     sprintf(levels,"%d",jac->levels);
290     args[cnt++] = (char*)"-level"; args[cnt++] = levels;
291   }
292   ierr = PetscOptionsTruth("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);CHKERRQ(ierr);
293   if (jac->bjilu) {
294     args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1";
295   }
296 
297   ierr = PetscOptionsTruth("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);CHKERRQ(ierr);
298   if (jac->printstatistics) {
299     args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1";
300     args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1";
301   }
302   ierr = PetscOptionsTail();CHKERRQ(ierr);
303   if (cnt) {
304     ierr = HYPRE_EuclidSetParams(jac->hsolver,cnt,args);CHKERRQ(ierr);
305   }
306   PetscFunctionReturn(0);
307 }
308 
309 #undef __FUNCT__
310 #define __FUNCT__ "PCView_HYPRE_Euclid"
311 static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
312 {
313   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
314   PetscErrorCode ierr;
315   PetscTruth     iascii;
316 
317   PetscFunctionBegin;
318   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
319   if (iascii) {
320     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n");CHKERRQ(ierr);
321     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr);
322     if (jac->bjilu) {
323       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr);
324     }
325   }
326   PetscFunctionReturn(0);
327 }
328 
329 /* --------------------------------------------------------------------------------------------*/
330 
331 #undef __FUNCT__
332 #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG"
333 static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
334 {
335   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
336   PetscErrorCode     ierr;
337   HYPRE_ParCSRMatrix hmat;
338   PetscScalar        *bv,*xv;
339   HYPRE_ParVector    jbv,jxv;
340   PetscScalar        *sbv,*sxv;
341   int                hierr;
342 
343   PetscFunctionBegin;
344   ierr = VecSet(x,0.0);CHKERRQ(ierr);
345   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
346   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
347   HYPREReplacePointer(jac->b,bv,sbv);
348   HYPREReplacePointer(jac->x,xv,sxv);
349 
350   ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr);
351   ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr);
352   ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr);
353 
354   hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv);
355   /* error code of 1 in BoomerAMG merely means convergence not achieved */
356   if (hierr && (hierr != 1)) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
357   if (hierr) hypre__global_error = 0;
358 
359   HYPREReplacePointer(jac->b,sbv,bv);
360   HYPREReplacePointer(jac->x,sxv,xv);
361   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
362   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
363   PetscFunctionReturn(0);
364 }
365 
366 static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
367 static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"};
368 static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
369 static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","","SOR/Jacobi","backward-SOR/Jacobi","","symmetric-SOR/Jacobi",
370                                                   "","","Gaussian-elimination"};
371 static const char *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i",
372                                                   "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"};
373 #undef __FUNCT__
374 #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG"
375 static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
376 {
377   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
378   PetscErrorCode ierr;
379   int            n,indx;
380   PetscTruth     flg, tmp_truth;
381   double         tmpdbl, twodbl[2];
382 
383   PetscFunctionBegin;
384   ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr);
385   ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr);
386   if (flg) {
387     jac->cycletype = indx;
388     ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr);
389   }
390   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr);
391   if (flg) {
392     if (jac->maxlevels < 2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
393     ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr);
394   }
395   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr);
396   if (flg) {
397     if (jac->maxiter < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
398     ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
399   }
400   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);
401   if (flg) {
402     if (jac->tol < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be greater than or equal to zero",jac->tol);
403     ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
404   }
405 
406   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr);
407   if (flg) {
408     if (jac->truncfactor < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor);
409     ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr);
410   }
411 
412  ierr = PetscOptionsInt("-pc_hypre_boomeramg_P_max","Max elements per row for interpolation operator ( 0=unlimited )","None",jac->pmax,&jac->pmax,&flg);CHKERRQ(ierr);
413   if (flg) {
414     if (jac->pmax < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"P_max %G must be greater than or equal to zero",jac->pmax);
415     ierr = HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);CHKERRQ(ierr);
416   }
417 
418  ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);CHKERRQ(ierr);
419   if (flg) {
420      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);
421 
422      ierr = HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl);CHKERRQ(ierr);
423   }
424 
425 
426  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);
427   if (flg) {
428      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);
429 
430      ierr = HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);CHKERRQ(ierr);
431   }
432 
433 
434   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr);
435   if (flg) {
436     if (jac->strongthreshold < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold);
437     ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr);
438   }
439   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr);
440   if (flg) {
441     if (jac->maxrowsum < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum);
442     if (jac->maxrowsum > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum);
443     ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr);
444   }
445 
446   /* Grid sweeps */
447   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);
448   if (flg) {
449     ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver,indx);CHKERRQ(ierr);
450     /* modify the jac structure so we can view the updated options with PC_View */
451     jac->gridsweeps[0] = indx;
452     jac->gridsweeps[1] = indx;
453     /*defaults coarse to 1 */
454     jac->gridsweeps[2] = 1;
455   }
456 
457   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx ,&flg);CHKERRQ(ierr);
458   if (flg) {
459     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 1);CHKERRQ(ierr);
460     jac->gridsweeps[0] = indx;
461   }
462   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr);
463   if (flg) {
464     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 2);CHKERRQ(ierr);
465     jac->gridsweeps[1] = indx;
466   }
467   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr);
468   if (flg) {
469     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 3);CHKERRQ(ierr);
470     jac->gridsweeps[2] = indx;
471   }
472 
473   /* Relax type */
474   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);
475   if (flg) {
476     jac->relaxtype[0] = jac->relaxtype[1]  = indx;
477     ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, indx);CHKERRQ(ierr);
478     /* by default, coarse type set to 9 */
479     jac->relaxtype[2] = 9;
480 
481   }
482   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
483   if (flg) {
484     jac->relaxtype[0] = indx;
485     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 1);CHKERRQ(ierr);
486   }
487   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
488   if (flg) {
489     jac->relaxtype[1] = indx;
490     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 2);CHKERRQ(ierr);
491   }
492   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr);
493   if (flg) {
494     jac->relaxtype[2] = indx;
495     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 3);CHKERRQ(ierr);
496   }
497 
498   /* Relaxation Weight */
499   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);
500   if (flg) {
501     ierr = HYPRE_BoomerAMGSetRelaxWt(jac->hsolver,tmpdbl);CHKERRQ(ierr);
502     jac->relaxweight = tmpdbl;
503   }
504 
505   n=2;
506   twodbl[0] = twodbl[1] = 1.0;
507   ierr = PetscOptionsRealArray("-pc_hypre_boomeramg_relax_weight_level","Set the relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);CHKERRQ(ierr);
508   if (flg) {
509     if (n == 2) {
510       indx =  (int)PetscAbsReal(twodbl[1]);
511       ierr = HYPRE_BoomerAMGSetLevelRelaxWt(jac->hsolver,twodbl[0],indx);CHKERRQ(ierr);
512     } else {
513       SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Relax weight level: you must provide 2 values separated by a comma (and no space), you provided %d",n);
514     }
515   }
516 
517   /* Outer relaxation Weight */
518   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);
519   if (flg) {
520     ierr = HYPRE_BoomerAMGSetOuterWt( jac->hsolver, tmpdbl);CHKERRQ(ierr);
521     jac->outerrelaxweight = tmpdbl;
522   }
523 
524   n=2;
525   twodbl[0] = twodbl[1] = 1.0;
526   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);
527   if (flg) {
528     if (n == 2) {
529       indx =  (int)PetscAbsReal(twodbl[1]);
530       ierr = HYPRE_BoomerAMGSetLevelOuterWt( jac->hsolver, twodbl[0], indx);CHKERRQ(ierr);
531     } else {
532       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);
533     }
534   }
535 
536   /* the Relax Order */
537   ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
538 
539   if (flg) {
540     jac->relaxorder = 0;
541     ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr);
542   }
543   ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr);
544   if (flg) {
545     jac->measuretype = indx;
546     ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr);
547   }
548   /* update list length 3/07 */
549   ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,11,HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr);
550   if (flg) {
551     jac->coarsentype = indx;
552     ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr);
553   }
554 
555   /* new 3/07 */
556   ierr = PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,14,HYPREBoomerAMGInterpType[0],&indx,&flg);CHKERRQ(ierr);
557   if (flg) {
558     jac->interptype = indx;
559     ierr = HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);CHKERRQ(ierr);
560   }
561 
562   flg  = PETSC_FALSE;
563   ierr = PetscOptionsTruth("-pc_hypre_boomeramg_print_statistics","Print statistics","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
564   if (flg) {
565     int level=3;
566     jac->printstatistics = PETSC_TRUE;
567     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
568     ierr = HYPRE_BoomerAMGSetPrintLevel(jac->hsolver,level);CHKERRQ(ierr);
569   }
570 
571   flg  = PETSC_FALSE;
572   ierr = PetscOptionsTruth("-pc_hypre_boomeramg_print_debug","Print debug information","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
573   if (flg) {
574     int level=3;
575     jac->printstatistics = PETSC_TRUE;
576     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
577     ierr = HYPRE_BoomerAMGSetDebugFlag(jac->hsolver,level);CHKERRQ(ierr);
578   }
579 
580   ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
581   if (flg && tmp_truth) {
582     jac->nodal_coarsen = 1;
583     ierr = HYPRE_BoomerAMGSetNodal(jac->hsolver,1);CHKERRQ(ierr);
584   }
585 
586   ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
587   if (flg && tmp_truth) {
588     PetscInt tmp_int;
589     ierr = PetscOptionsInt( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);CHKERRQ(ierr);
590     if (flg) jac->nodal_relax_levels = tmp_int;
591     ierr = HYPRE_BoomerAMGSetSmoothType(jac->hsolver,6);CHKERRQ(ierr);
592     ierr = HYPRE_BoomerAMGSetDomainType(jac->hsolver,1);CHKERRQ(ierr);
593     ierr = HYPRE_BoomerAMGSetOverlap(jac->hsolver,0);CHKERRQ(ierr);
594     ierr = HYPRE_BoomerAMGSetSmoothNumLevels(jac->hsolver,jac->nodal_relax_levels);CHKERRQ(ierr);
595   }
596 
597   ierr = PetscOptionsTail();CHKERRQ(ierr);
598   PetscFunctionReturn(0);
599 }
600 
601 #undef __FUNCT__
602 #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG"
603 static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
604 {
605   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
606   PetscErrorCode ierr;
607   int            oits;
608 
609   PetscFunctionBegin;
610   ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,its*jac->maxiter);CHKERRQ(ierr);
611   ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,rtol);CHKERRQ(ierr);
612   jac->applyrichardson = PETSC_TRUE;
613   ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr);
614   jac->applyrichardson = PETSC_FALSE;
615   ierr = HYPRE_BoomerAMGGetNumIterations(jac->hsolver,&oits);CHKERRQ(ierr);
616   *outits = oits;
617   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
618   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
619   ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
620   ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
621   PetscFunctionReturn(0);
622 }
623 
624 
625 #undef __FUNCT__
626 #define __FUNCT__ "PCView_HYPRE_BoomerAMG"
627 static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
628 {
629   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
630   PetscErrorCode ierr;
631   PetscTruth     iascii;
632 
633   PetscFunctionBegin;
634   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
635   if (iascii) {
636     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr);
637     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr);
638     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr);
639     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr);
640     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr);
641     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr);
642     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation truncation factor %G\n",jac->truncfactor);CHKERRQ(ierr);
643     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);CHKERRQ(ierr);
644     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);CHKERRQ(ierr);
645     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);CHKERRQ(ierr);
646 
647     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr);
648 
649     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps down         %d\n",jac->gridsweeps[0]);CHKERRQ(ierr);
650     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps up           %d\n",jac->gridsweeps[1]);CHKERRQ(ierr);
651     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on coarse    %d\n",jac->gridsweeps[2]);CHKERRQ(ierr);
652 
653     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr);
654     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr);
655     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr);
656 
657     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax weight  (all)      %G\n",jac->relaxweight);CHKERRQ(ierr);
658     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr);
659 
660     if (jac->relaxorder) {
661       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr);
662     } else {
663       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr);
664     }
665     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr);
666     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr);
667     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation type  %s\n",HYPREBoomerAMGInterpType[jac->interptype]);CHKERRQ(ierr);
668     if (jac->nodal_coarsen) {
669       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal())\n");CHKERRQ(ierr);
670     }
671     if (jac->nodal_relax) {
672       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);CHKERRQ(ierr);
673     }
674   }
675   PetscFunctionReturn(0);
676 }
677 
678 /* --------------------------------------------------------------------------------------------*/
679 #undef __FUNCT__
680 #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails"
681 static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc)
682 {
683   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
684   PetscErrorCode ierr;
685   int            indx;
686   PetscTruth     flag;
687   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};
688 
689   PetscFunctionBegin;
690   ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr);
691   ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr);
692   ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr);
693   if (flag) {
694     ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr);
695   }
696 
697   ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr);
698   if (flag) {
699     ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr);
700   }
701 
702   ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr);
703   if (flag) {
704     ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr);
705   }
706 
707   ierr = PetscOptionsTruth("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,&flag);CHKERRQ(ierr);
708   if (flag) {
709     ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr);
710   }
711 
712   ierr = PetscOptionsTruth("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,&flag);CHKERRQ(ierr);
713   if (flag) {
714     ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr);
715   }
716 
717   ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);CHKERRQ(ierr);
718   if (flag) {
719     jac->symt = indx;
720     ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr);
721   }
722 
723   ierr = PetscOptionsTail();CHKERRQ(ierr);
724   PetscFunctionReturn(0);
725 }
726 
727 #undef __FUNCT__
728 #define __FUNCT__ "PCView_HYPRE_ParaSails"
729 static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
730 {
731   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
732   PetscErrorCode ierr;
733   PetscTruth     iascii;
734   const char     *symt = 0;;
735 
736   PetscFunctionBegin;
737   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
738   if (iascii) {
739     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");CHKERRQ(ierr);
740     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr);
741     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr);
742     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr);
743     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr);
744     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscTruths[jac->ruse]);CHKERRQ(ierr);
745     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscTruths[jac->logging]);CHKERRQ(ierr);
746     if (!jac->symt) {
747       symt = "nonsymmetric matrix and preconditioner";
748     } else if (jac->symt == 1) {
749       symt = "SPD matrix and preconditioner";
750     } else if (jac->symt == 2) {
751       symt = "nonsymmetric matrix but SPD preconditioner";
752     } else {
753       SETERRQ1(PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
754     }
755     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr);
756   }
757   PetscFunctionReturn(0);
758 }
759 /* ---------------------------------------------------------------------------------*/
760 
761 EXTERN_C_BEGIN
762 #undef __FUNCT__
763 #define __FUNCT__ "PCHYPREGetType_HYPRE"
764 PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType_HYPRE(PC pc,const char *name[])
765 {
766   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
767 
768   PetscFunctionBegin;
769   *name = jac->hypre_type;
770   PetscFunctionReturn(0);
771 }
772 EXTERN_C_END
773 
774 EXTERN_C_BEGIN
775 #undef __FUNCT__
776 #define __FUNCT__ "PCHYPRESetType_HYPRE"
777 PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType_HYPRE(PC pc,const char name[])
778 {
779   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
780   PetscErrorCode ierr;
781   PetscTruth     flag;
782 
783   PetscFunctionBegin;
784   if (jac->hypre_type) {
785     ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr);
786     if (!flag) SETERRQ(PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
787     PetscFunctionReturn(0);
788   } else {
789     ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr);
790   }
791 
792   jac->maxiter            = PETSC_DEFAULT;
793   jac->tol                = PETSC_DEFAULT;
794   jac->printstatistics    = PetscLogPrintInfo;
795 
796   ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr);
797   if (flag) {
798     ierr                    = HYPRE_ParCSRPilutCreate(jac->comm_hypre,&jac->hsolver);
799     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
800     pc->ops->view           = PCView_HYPRE_Pilut;
801     jac->destroy            = HYPRE_ParCSRPilutDestroy;
802     jac->setup              = HYPRE_ParCSRPilutSetup;
803     jac->solve              = HYPRE_ParCSRPilutSolve;
804     jac->factorrowsize      = PETSC_DEFAULT;
805     PetscFunctionReturn(0);
806   }
807   ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr);
808   if (flag) {
809     ierr                    = HYPRE_ParaSailsCreate(jac->comm_hypre,&jac->hsolver);
810     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
811     pc->ops->view           = PCView_HYPRE_ParaSails;
812     jac->destroy            = HYPRE_ParaSailsDestroy;
813     jac->setup              = HYPRE_ParaSailsSetup;
814     jac->solve              = HYPRE_ParaSailsSolve;
815     /* initialize */
816     jac->nlevels     = 1;
817     jac->threshhold  = .1;
818     jac->filter      = .1;
819     jac->loadbal     = 0;
820     if (PetscLogPrintInfo) {
821       jac->logging     = (int) PETSC_TRUE;
822     } else {
823       jac->logging     = (int) PETSC_FALSE;
824     }
825     jac->ruse = (int) PETSC_FALSE;
826     jac->symt   = 0;
827     ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr);
828     ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr);
829     ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr);
830     ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr);
831     ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr);
832     ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr);
833     PetscFunctionReturn(0);
834   }
835   ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr);
836   if (flag) {
837     ierr                    = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver);
838     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
839     pc->ops->view           = PCView_HYPRE_Euclid;
840     jac->destroy            = HYPRE_EuclidDestroy;
841     jac->setup              = HYPRE_EuclidSetup;
842     jac->solve              = HYPRE_EuclidSolve;
843     /* initialization */
844     jac->bjilu              = PETSC_FALSE;
845     jac->levels             = 1;
846     PetscFunctionReturn(0);
847   }
848   ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr);
849   if (flag) {
850     ierr                     = HYPRE_BoomerAMGCreate(&jac->hsolver);
851     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
852     pc->ops->view            = PCView_HYPRE_BoomerAMG;
853     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
854     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
855     jac->destroy             = HYPRE_BoomerAMGDestroy;
856     jac->setup               = HYPRE_BoomerAMGSetup;
857     jac->solve               = HYPRE_BoomerAMGSolve;
858     jac->applyrichardson     = PETSC_FALSE;
859     /* these defaults match the hypre defaults */
860     jac->cycletype        = 1;
861     jac->maxlevels        = 25;
862     jac->maxiter          = 1;
863     jac->tol              = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */
864     jac->truncfactor      = 0.0;
865     jac->strongthreshold  = .25;
866     jac->maxrowsum        = .9;
867     jac->coarsentype      = 6;
868     jac->measuretype      = 0;
869     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = 1;
870     jac->relaxtype[0]     = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */
871     jac->relaxtype[2]     = 9; /*G.E. */
872     jac->relaxweight      = 1.0;
873     jac->outerrelaxweight = 1.0;
874     jac->relaxorder       = 1;
875     jac->interptype       = 0;
876     jac->agg_nl           = 0;
877     jac->pmax             = 0;
878     jac->truncfactor      = 0.0;
879     jac->agg_num_paths    = 1;
880 
881     jac->nodal_coarsen    = 0;
882     jac->nodal_relax      = PETSC_FALSE;
883     jac->nodal_relax_levels = 1;
884     ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr);
885     ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr);
886     ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
887     ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
888     ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr);
889     ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr);
890     ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr);
891     ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr);
892     ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr);
893     ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr);
894     ierr = HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);CHKERRQ(ierr);
895     ierr = HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl);
896     ierr = HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);CHKERRQ(ierr);
897     ierr = HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);CHKERRQ(ierr);
898     ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, jac->relaxtype[0]);  /*defaults coarse to 9*/
899     ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver, jac->gridsweeps[0]); /*defaults coarse to 1 */
900     PetscFunctionReturn(0);
901   }
902   ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr);
903   jac->hypre_type = PETSC_NULL;
904   SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
905   PetscFunctionReturn(0);
906 }
907 EXTERN_C_END
908 
909 /*
910     It only gets here if the HYPRE type has not been set before the call to
911    ...SetFromOptions() which actually is most of the time
912 */
913 #undef __FUNCT__
914 #define __FUNCT__ "PCSetFromOptions_HYPRE"
915 static PetscErrorCode PCSetFromOptions_HYPRE(PC pc)
916 {
917   PetscErrorCode ierr;
918   int            indx;
919   const char     *type[] = {"pilut","parasails","boomeramg","euclid"};
920   PetscTruth     flg;
921 
922   PetscFunctionBegin;
923   ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr);
924   ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"boomeramg",&indx,&flg);CHKERRQ(ierr);
925   if (flg) {
926     ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr);
927   } else {
928     ierr = PCHYPRESetType_HYPRE(pc,"boomeramg");CHKERRQ(ierr);
929   }
930   if (pc->ops->setfromoptions) {
931     ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr);
932   }
933   ierr = PetscOptionsTail();CHKERRQ(ierr);
934   PetscFunctionReturn(0);
935 }
936 
937 #undef __FUNCT__
938 #define __FUNCT__ "PCHYPRESetType"
939 /*@C
940      PCHYPRESetType - Sets which hypre preconditioner you wish to use
941 
942    Input Parameters:
943 +     pc - the preconditioner context
944 -     name - either  pilut, parasails, boomeramg, euclid
945 
946    Options Database Keys:
947    -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
948 
949    Level: intermediate
950 
951 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
952            PCHYPRE
953 
954 @*/
955 PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC pc,const char name[])
956 {
957   PetscErrorCode ierr,(*f)(PC,const char[]);
958 
959   PetscFunctionBegin;
960   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
961   PetscValidCharPointer(name,2);
962   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPRESetType_C",(void (**)(void))&f);CHKERRQ(ierr);
963   if (f) {
964     ierr = (*f)(pc,name);CHKERRQ(ierr);
965   }
966   PetscFunctionReturn(0);
967 }
968 
969 #undef __FUNCT__
970 #define __FUNCT__ "PCHYPREGetType"
971 /*@C
972      PCHYPREGetType - Gets which hypre preconditioner you are using
973 
974    Input Parameter:
975 .     pc - the preconditioner context
976 
977    Output Parameter:
978 .     name - either  pilut, parasails, boomeramg, euclid
979 
980    Level: intermediate
981 
982 .seealso:  PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC,
983            PCHYPRE
984 
985 @*/
986 PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC pc,const char *name[])
987 {
988   PetscErrorCode ierr,(*f)(PC,const char*[]);
989 
990   PetscFunctionBegin;
991   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
992   PetscValidPointer(name,2);
993   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPREGetType_C",(void (**)(void))&f);CHKERRQ(ierr);
994   if (f) {
995     ierr = (*f)(pc,name);CHKERRQ(ierr);
996   }
997   PetscFunctionReturn(0);
998 }
999 
1000 /*MC
1001      PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre
1002 
1003    Options Database Keys:
1004 +   -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
1005 -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
1006           preconditioner
1007 
1008    Level: intermediate
1009 
1010    Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()),
1011           the many hypre options can ONLY be set via the options database (e.g. the command line
1012           or with PetscOptionsSetValue(), there are no functions to set them)
1013 
1014           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations
1015           (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if
1016           -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner
1017           (-pc_hypre_boomeramg_rtol should be set to 0.0 - the default - to strictly use a fixed number of
1018           iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations
1019           and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10
1020           then AT MOST twenty V-cycles of boomeramg will be called.
1021 
1022            Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation
1023            (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry.
1024            Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi.
1025           If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly
1026           and use -ksp_max_it to control the number of V-cycles.
1027           (see the PETSc FAQ.html at the PETSc website under the Documentation tab).
1028 
1029 	  2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option
1030 	  -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L.
1031 
1032           See PCPFMG for access to the hypre Struct PFMG solver
1033 
1034 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1035            PCHYPRESetType(), PCPFMG
1036 
1037 M*/
1038 
1039 EXTERN_C_BEGIN
1040 #undef __FUNCT__
1041 #define __FUNCT__ "PCCreate_HYPRE"
1042 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_HYPRE(PC pc)
1043 {
1044   PC_HYPRE       *jac;
1045   PetscErrorCode ierr;
1046 
1047   PetscFunctionBegin;
1048   ierr = PetscNewLog(pc,PC_HYPRE,&jac);CHKERRQ(ierr);
1049   pc->data                 = jac;
1050   pc->ops->destroy         = PCDestroy_HYPRE;
1051   pc->ops->setfromoptions  = PCSetFromOptions_HYPRE;
1052   pc->ops->setup           = PCSetUp_HYPRE;
1053   pc->ops->apply           = PCApply_HYPRE;
1054   jac->comm_hypre          = MPI_COMM_NULL;
1055   jac->hypre_type          = PETSC_NULL;
1056   /* duplicate communicator for hypre */
1057   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(jac->comm_hypre));CHKERRQ(ierr);
1058   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr);
1059   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr);
1060   PetscFunctionReturn(0);
1061 }
1062 EXTERN_C_END
1063 
1064 /* ---------------------------------------------------------------------------------------------------------------------------------*/
1065 
1066 /* working with a HYPRE_StructMatrix and need access to its data */
1067 #include "../src/dm/da/utils/mhyp.h"
1068 /* this include is needed ONLY to allow access to the private data inside the Mat object specific to hypre */
1069 #include "private/matimpl.h"
1070 
1071 typedef struct {
1072   MPI_Comm            hcomm;       /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */
1073   HYPRE_StructSolver  hsolver;
1074 
1075   /* keep copy of PFMG options used so may view them */
1076   int                 its;
1077   double              tol;
1078   int                 relax_type;
1079   int                 rap_type;
1080   int                 num_pre_relax,num_post_relax;
1081   int                 max_levels;
1082 } PC_PFMG;
1083 
1084 #undef __FUNCT__
1085 #define __FUNCT__ "PCDestroy_PFMG"
1086 PetscErrorCode PCDestroy_PFMG(PC pc)
1087 {
1088   PetscErrorCode ierr;
1089   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1090 
1091   PetscFunctionBegin;
1092   if (ex->hsolver) {ierr = HYPRE_StructPFMGDestroy(ex->hsolver);CHKERRQ(ierr);}
1093   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1094   ierr = PetscFree(ex);CHKERRQ(ierr);
1095   PetscFunctionReturn(0);
1096 }
1097 
1098 static const char *PFMGRelaxType[]   = {"Jacobi","Weighted-Jacobi","symmetric-Red/Black-Gauss-Seidel","Red/Black-Gauss-Seidel"};
1099 static const char *PFMGRAPType[]   = {"Galerkin","non-Galerkin"};
1100 
1101 #undef __FUNCT__
1102 #define __FUNCT__ "PCView_PFMG"
1103 PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer)
1104 {
1105   PetscErrorCode ierr;
1106   PetscTruth     iascii;
1107   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1108 
1109   PetscFunctionBegin;
1110   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1111   if (iascii) {
1112     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG preconditioning\n");CHKERRQ(ierr);
1113     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
1114     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
1115     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
1116     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: RAP type %s\n",PFMGRAPType[ex->rap_type]);CHKERRQ(ierr);
1117     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
1118     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max levels %d\n",ex->max_levels);CHKERRQ(ierr);
1119   }
1120   PetscFunctionReturn(0);
1121 }
1122 
1123 
1124 #undef __FUNCT__
1125 #define __FUNCT__ "PCSetFromOptions_PFMG"
1126 PetscErrorCode PCSetFromOptions_PFMG(PC pc)
1127 {
1128   PetscErrorCode ierr;
1129   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1130   PetscTruth     flg = PETSC_FALSE;
1131 
1132   PetscFunctionBegin;
1133   ierr = PetscOptionsHead("PFMG options");CHKERRQ(ierr);
1134   ierr = PetscOptionsTruth("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1135   if (flg) {
1136     int level=3;
1137     ierr = HYPRE_StructPFMGSetPrintLevel(ex->hsolver,level);CHKERRQ(ierr);
1138   }
1139   ierr = PetscOptionsInt("-pc_pfmg_its","Number of iterations of PFMG to use as preconditioner","HYPRE_StructPFMGSetMaxIter",ex->its,&ex->its,PETSC_NULL);CHKERRQ(ierr);
1140   ierr = HYPRE_StructPFMGSetMaxIter(ex->hsolver,ex->its);CHKERRQ(ierr);
1141   ierr = PetscOptionsInt("-pc_pfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_StructPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,PETSC_NULL);CHKERRQ(ierr);
1142   ierr = HYPRE_StructPFMGSetNumPreRelax(ex->hsolver,ex->num_pre_relax);CHKERRQ(ierr);
1143   ierr = PetscOptionsInt("-pc_pfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_StructPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,PETSC_NULL);CHKERRQ(ierr);
1144   ierr = HYPRE_StructPFMGSetNumPostRelax(ex->hsolver,ex->num_post_relax);CHKERRQ(ierr);
1145 
1146   ierr = PetscOptionsInt("-pc_pfmg_max_levels","Max Levels for MG hierarchy","HYPRE_StructPFMGSetMaxLevels",ex->max_levels,&ex->max_levels,PETSC_NULL);CHKERRQ(ierr);
1147   ierr = HYPRE_StructPFMGSetMaxLevels(ex->hsolver,ex->max_levels);CHKERRQ(ierr);
1148 
1149   ierr = PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr);
1150   ierr = HYPRE_StructPFMGSetTol(ex->hsolver,ex->tol);CHKERRQ(ierr);
1151   ierr = PetscOptionsEList("-pc_pfmg_relax_type","Relax type for the up and down cycles","HYPRE_StructPFMGSetRelaxType",PFMGRelaxType,4,PFMGRelaxType[ex->relax_type],&ex->relax_type,PETSC_NULL);CHKERRQ(ierr);
1152   ierr = HYPRE_StructPFMGSetRelaxType(ex->hsolver, ex->relax_type);CHKERRQ(ierr);
1153   ierr = PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,2,PFMGRAPType[ex->rap_type],&ex->rap_type,PETSC_NULL);CHKERRQ(ierr);
1154   ierr = HYPRE_StructPFMGSetRAPType(ex->hsolver, ex->rap_type);CHKERRQ(ierr);
1155   ierr = PetscOptionsTail();CHKERRQ(ierr);
1156   PetscFunctionReturn(0);
1157 }
1158 
1159 #undef __FUNCT__
1160 #define __FUNCT__ "PCApply_PFMG"
1161 PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y)
1162 {
1163   PetscErrorCode  ierr;
1164   PC_PFMG         *ex = (PC_PFMG*) pc->data;
1165   PetscScalar     *xx,*yy;
1166   int             ilower[3],iupper[3];
1167   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
1168 
1169   PetscFunctionBegin;
1170   ierr = DAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1171   iupper[0] += ilower[0] - 1;
1172   iupper[1] += ilower[1] - 1;
1173   iupper[2] += ilower[2] - 1;
1174 
1175   /* copy x values over to hypre */
1176   ierr = HYPRE_StructVectorSetConstantValues(mx->hb,0.0);CHKERRQ(ierr);
1177   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1178   ierr = HYPRE_StructVectorSetBoxValues(mx->hb,ilower,iupper,xx);CHKERRQ(ierr);
1179   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1180   ierr = HYPRE_StructVectorAssemble(mx->hb);CHKERRQ(ierr);
1181 
1182   ierr = HYPRE_StructPFMGSolve(ex->hsolver,mx->hmat,mx->hb,mx->hx);CHKERRQ(ierr);
1183 
1184   /* copy solution values back to PETSc */
1185   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1186   ierr = HYPRE_StructVectorGetBoxValues(mx->hx,ilower,iupper,yy);CHKERRQ(ierr);
1187   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1188   PetscFunctionReturn(0);
1189 }
1190 
1191 #undef __FUNCT__
1192 #define __FUNCT__ "PCApplyRichardson_PFMG"
1193 static PetscErrorCode PCApplyRichardson_PFMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
1194 {
1195   PC_PFMG        *jac = (PC_PFMG*)pc->data;
1196   PetscErrorCode ierr;
1197   int            oits;
1198 
1199   PetscFunctionBegin;
1200   ierr = HYPRE_StructPFMGSetMaxIter(jac->hsolver,its*jac->its);CHKERRQ(ierr);
1201   ierr = HYPRE_StructPFMGSetTol(jac->hsolver,rtol);CHKERRQ(ierr);
1202 
1203   ierr = PCApply_PFMG(pc,b,y);CHKERRQ(ierr);
1204   ierr = HYPRE_StructPFMGGetNumIterations(jac->hsolver,&oits);CHKERRQ(ierr);
1205   *outits = oits;
1206   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
1207   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
1208   ierr = HYPRE_StructPFMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
1209   ierr = HYPRE_StructPFMGSetMaxIter(jac->hsolver,jac->its);CHKERRQ(ierr);
1210   PetscFunctionReturn(0);
1211 }
1212 
1213 
1214 #undef __FUNCT__
1215 #define __FUNCT__ "PCSetUp_PFMG"
1216 PetscErrorCode PCSetUp_PFMG(PC pc)
1217 {
1218   PetscErrorCode  ierr;
1219   PC_PFMG         *ex = (PC_PFMG*) pc->data;
1220   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
1221   PetscTruth      flg;
1222 
1223   PetscFunctionBegin;
1224   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATHYPRESTRUCT,&flg);CHKERRQ(ierr);
1225   if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESTRUCT with this preconditioner");
1226 
1227   /* create the hypre solver object and set its information */
1228   if (ex->hsolver) {
1229     ierr = HYPRE_StructPFMGDestroy(ex->hsolver);CHKERRQ(ierr);
1230   }
1231   ierr = HYPRE_StructPFMGCreate(ex->hcomm,&ex->hsolver);CHKERRQ(ierr);
1232   ierr = PCSetFromOptions_PFMG(pc);CHKERRQ(ierr);
1233   ierr = HYPRE_StructPFMGSetup(ex->hsolver,mx->hmat,mx->hb,mx->hx);CHKERRQ(ierr);
1234   ierr = HYPRE_StructPFMGSetZeroGuess(ex->hsolver);CHKERRQ(ierr);
1235 
1236   PetscFunctionReturn(0);
1237 }
1238 
1239 
1240 /*MC
1241      PCPFMG - the hypre PFMG multigrid solver
1242 
1243    Level: advanced
1244 
1245    Options Database:
1246 + -pc_pfmg_its <its> number of iterations of PFMG to use as preconditioner
1247 . -pc_pfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
1248 . -pc_pfmg_num_post_relax <steps> number of smoothing steps after coarse grid
1249 . -pc_pfmg_tol <tol> tolerance of PFMG
1250 . -pc_pfmg_relax_type -relaxation type for the up and down cycles, one of Jacobi,Weighted-Jacobi,symmetric-Red/Black-Gauss-Seidel,Red/Black-Gauss-Seidel
1251 - -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin
1252 
1253    Notes:  This is for CELL-centered descretizations
1254 
1255            This must be used with the MATHYPRESTRUCT matrix type.
1256            This is less general than in hypre, it supports only one block per process defined by a PETSc DA.
1257 
1258 .seealso:  PCMG, MATHYPRESTRUCT
1259 M*/
1260 
1261 EXTERN_C_BEGIN
1262 #undef __FUNCT__
1263 #define __FUNCT__ "PCCreate_PFMG"
1264 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_PFMG(PC pc)
1265 {
1266   PetscErrorCode ierr;
1267   PC_PFMG        *ex;
1268 
1269   PetscFunctionBegin;
1270   ierr = PetscNew(PC_PFMG,&ex);CHKERRQ(ierr);\
1271   pc->data = ex;
1272 
1273   ex->its            = 1;
1274   ex->tol            = 1.e-8;
1275   ex->relax_type     = 1;
1276   ex->rap_type       = 0;
1277   ex->num_pre_relax  = 1;
1278   ex->num_post_relax = 1;
1279   ex->max_levels     = 0;
1280 
1281   pc->ops->setfromoptions  = PCSetFromOptions_PFMG;
1282   pc->ops->view            = PCView_PFMG;
1283   pc->ops->destroy         = PCDestroy_PFMG;
1284   pc->ops->apply           = PCApply_PFMG;
1285   pc->ops->applyrichardson = PCApplyRichardson_PFMG;
1286   pc->ops->setup           = PCSetUp_PFMG;
1287   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr);
1288   ierr = HYPRE_StructPFMGCreate(ex->hcomm,&ex->hsolver);CHKERRQ(ierr);
1289   PetscFunctionReturn(0);
1290 }
1291 EXTERN_C_END
1292 
1293 /* we know we are working with a HYPRE_SStructMatrix */
1294 typedef struct {
1295   MPI_Comm            hcomm;       /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */
1296   HYPRE_SStructSolver  ss_solver;
1297 
1298   /* keep copy of SYSPFMG options used so may view them */
1299   int                 its;
1300   double              tol;
1301   int                 relax_type;
1302   int                 num_pre_relax,num_post_relax;
1303 } PC_SysPFMG;
1304 
1305 #undef __FUNCT__
1306 #define __FUNCT__ "PCDestroy_SysPFMG"
1307 PetscErrorCode PCDestroy_SysPFMG(PC pc)
1308 {
1309   PetscErrorCode ierr;
1310   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1311 
1312   PetscFunctionBegin;
1313   if (ex->ss_solver) {ierr = HYPRE_SStructSysPFMGDestroy(ex->ss_solver);CHKERRQ(ierr);}
1314   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1315   ierr = PetscFree(ex);CHKERRQ(ierr);
1316   PetscFunctionReturn(0);
1317 }
1318 
1319 static const char *SysPFMGRelaxType[]   = {"Weighted-Jacobi","Red/Black-Gauss-Seidel"};
1320 
1321 #undef __FUNCT__
1322 #define __FUNCT__ "PCView_SysPFMG"
1323 PetscErrorCode PCView_SysPFMG(PC pc,PetscViewer viewer)
1324 {
1325   PetscErrorCode ierr;
1326   PetscTruth     iascii;
1327   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1328 
1329   PetscFunctionBegin;
1330   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
1331   if (iascii) {
1332     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG preconditioning\n");CHKERRQ(ierr);
1333     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
1334     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
1335     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
1336     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
1337   }
1338   PetscFunctionReturn(0);
1339 }
1340 
1341 
1342 #undef __FUNCT__
1343 #define __FUNCT__ "PCSetFromOptions_SysPFMG"
1344 PetscErrorCode PCSetFromOptions_SysPFMG(PC pc)
1345 {
1346   PetscErrorCode ierr;
1347   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1348   PetscTruth     flg = PETSC_FALSE;
1349 
1350   PetscFunctionBegin;
1351   ierr = PetscOptionsHead("SysPFMG options");CHKERRQ(ierr);
1352   ierr = PetscOptionsTruth("-pc_syspfmg_print_statistics","Print statistics","HYPRE_SStructSysPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1353   if (flg) {
1354     int level=3;
1355     ierr = HYPRE_SStructSysPFMGSetPrintLevel(ex->ss_solver,level);CHKERRQ(ierr);
1356   }
1357   ierr = PetscOptionsInt("-pc_syspfmg_its","Number of iterations of SysPFMG to use as preconditioner","HYPRE_SStructSysPFMGSetMaxIter",ex->its,&ex->its,PETSC_NULL);CHKERRQ(ierr);
1358   ierr = HYPRE_SStructSysPFMGSetMaxIter(ex->ss_solver,ex->its);CHKERRQ(ierr);
1359   ierr = PetscOptionsInt("-pc_syspfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_SStructSysPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,PETSC_NULL);CHKERRQ(ierr);
1360   ierr = HYPRE_SStructSysPFMGSetNumPreRelax(ex->ss_solver,ex->num_pre_relax);CHKERRQ(ierr);
1361   ierr = PetscOptionsInt("-pc_syspfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_SStructSysPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,PETSC_NULL);CHKERRQ(ierr);
1362   ierr = HYPRE_SStructSysPFMGSetNumPostRelax(ex->ss_solver,ex->num_post_relax);CHKERRQ(ierr);
1363 
1364   ierr = PetscOptionsReal("-pc_syspfmg_tol","Tolerance of SysPFMG","HYPRE_SStructSysPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr);
1365   ierr = HYPRE_SStructSysPFMGSetTol(ex->ss_solver,ex->tol);CHKERRQ(ierr);
1366   ierr = PetscOptionsEList("-pc_syspfmg_relax_type","Relax type for the up and down cycles","HYPRE_SStructSysPFMGSetRelaxType",SysPFMGRelaxType,4,SysPFMGRelaxType[ex->relax_type],&ex->relax_type,PETSC_NULL);CHKERRQ(ierr);
1367   ierr = HYPRE_SStructSysPFMGSetRelaxType(ex->ss_solver, ex->relax_type);CHKERRQ(ierr);
1368   ierr = PetscOptionsTail();CHKERRQ(ierr);
1369   PetscFunctionReturn(0);
1370 }
1371 
1372 #undef __FUNCT__
1373 #define __FUNCT__ "PCApply_SysPFMG"
1374 PetscErrorCode PCApply_SysPFMG(PC pc,Vec x,Vec y)
1375 {
1376   PetscErrorCode    ierr;
1377   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
1378   PetscScalar      *xx,*yy;
1379   int               ilower[3],iupper[3];
1380   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct *)(pc->pmat->data);
1381   int               ordering= mx->dofs_order;
1382   int               nvars= mx->nvars;
1383   int               part= 0;
1384   int               size;
1385   int               i;
1386 
1387   PetscFunctionBegin;
1388   ierr = DAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1389   iupper[0] += ilower[0] - 1;
1390   iupper[1] += ilower[1] - 1;
1391   iupper[2] += ilower[2] - 1;
1392 
1393   size= 1;
1394   for (i= 0; i< 3; i++) {
1395      size*= (iupper[i]-ilower[i]+1);
1396   }
1397   /* copy x values over to hypre for variable ordering */
1398   if (ordering) {
1399      ierr = HYPRE_SStructVectorSetConstantValues(mx->ss_b,0.0);CHKERRQ(ierr);
1400      ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1401      for (i= 0; i< nvars; i++) {
1402         ierr = HYPRE_SStructVectorSetBoxValues(mx->ss_b,part,ilower,iupper,i,xx+(size*i));CHKERRQ(ierr);
1403      }
1404      ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1405      ierr = HYPRE_SStructVectorAssemble(mx->ss_b);CHKERRQ(ierr);
1406 
1407      ierr = HYPRE_SStructMatrixMatvec(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x);CHKERRQ(ierr);
1408      ierr = HYPRE_SStructSysPFMGSolve(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x);CHKERRQ(ierr);
1409 
1410      /* copy solution values back to PETSc */
1411      ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1412      for (i= 0; i< nvars; i++) {
1413         ierr = HYPRE_SStructVectorGetBoxValues(mx->ss_x,part,ilower,iupper,i,yy+(size*i));CHKERRQ(ierr);
1414      }
1415      ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1416   }
1417 
1418   else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
1419      PetscScalar     *z;
1420      int              j, k;
1421 
1422      ierr = PetscMalloc(nvars*size*sizeof(PetscScalar),&z);CHKERRQ(ierr);
1423      ierr = HYPRE_SStructVectorSetConstantValues(mx->ss_b,0.0);CHKERRQ(ierr);
1424      ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1425 
1426      /* transform nodal to hypre's variable ordering for sys_pfmg */
1427      for (i= 0; i< size; i++) {
1428         k= i*nvars;
1429         for (j= 0; j< nvars; j++) {
1430            z[j*size+i]= xx[k+j];
1431         }
1432      }
1433      for (i= 0; i< nvars; i++) {
1434         ierr = HYPRE_SStructVectorSetBoxValues(mx->ss_b,part,ilower,iupper,i,z+(size*i));CHKERRQ(ierr);
1435      }
1436      ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1437 
1438      ierr = HYPRE_SStructVectorAssemble(mx->ss_b);CHKERRQ(ierr);
1439 
1440      ierr = HYPRE_SStructSysPFMGSolve(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x);CHKERRQ(ierr);
1441 
1442      /* copy solution values back to PETSc */
1443      ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1444      for (i= 0; i< nvars; i++) {
1445         ierr = HYPRE_SStructVectorGetBoxValues(mx->ss_x,part,ilower,iupper,i,z+(size*i));CHKERRQ(ierr);
1446      }
1447      /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
1448      for (i= 0; i< size; i++) {
1449         k= i*nvars;
1450         for (j= 0; j< nvars; j++) {
1451            yy[k+j]= z[j*size+i];
1452         }
1453      }
1454      ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1455 
1456      ierr = PetscFree(z);CHKERRQ(ierr);
1457   }
1458 
1459   PetscFunctionReturn(0);
1460 }
1461 
1462 #undef __FUNCT__
1463 #define __FUNCT__ "PCApplyRichardson_SysPFMG"
1464 static PetscErrorCode PCApplyRichardson_SysPFMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
1465 {
1466   PC_SysPFMG    *jac = (PC_SysPFMG*)pc->data;
1467   PetscErrorCode ierr;
1468   int            oits;
1469 
1470   PetscFunctionBegin;
1471 
1472   ierr = HYPRE_SStructSysPFMGSetMaxIter(jac->ss_solver,its*jac->its);CHKERRQ(ierr);
1473   ierr = HYPRE_SStructSysPFMGSetTol(jac->ss_solver,rtol);CHKERRQ(ierr);
1474 
1475   ierr = PCApply_SysPFMG(pc,b,y);CHKERRQ(ierr);
1476   ierr = HYPRE_SStructSysPFMGGetNumIterations(jac->ss_solver,&oits);CHKERRQ(ierr);
1477   *outits = oits;
1478   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
1479   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
1480   ierr = HYPRE_SStructSysPFMGSetTol(jac->ss_solver,jac->tol);CHKERRQ(ierr);
1481   ierr = HYPRE_SStructSysPFMGSetMaxIter(jac->ss_solver,jac->its);CHKERRQ(ierr);
1482   PetscFunctionReturn(0);
1483 }
1484 
1485 
1486 #undef __FUNCT__
1487 #define __FUNCT__ "PCSetUp_SysPFMG"
1488 PetscErrorCode PCSetUp_SysPFMG(PC pc)
1489 {
1490   PetscErrorCode    ierr;
1491   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
1492   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct *)(pc->pmat->data);
1493   PetscTruth        flg;
1494 
1495   PetscFunctionBegin;
1496   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATHYPRESSTRUCT,&flg);CHKERRQ(ierr);
1497   if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESSTRUCT with this preconditioner");
1498 
1499   /* create the hypre sstruct solver object and set its information */
1500   if (ex->ss_solver) {
1501     ierr = HYPRE_SStructSysPFMGDestroy(ex->ss_solver);CHKERRQ(ierr);
1502   }
1503   ierr = HYPRE_SStructSysPFMGCreate(ex->hcomm,&ex->ss_solver);CHKERRQ(ierr);
1504   ierr = PCSetFromOptions_SysPFMG(pc);CHKERRQ(ierr);
1505   ierr = HYPRE_SStructSysPFMGSetZeroGuess(ex->ss_solver);CHKERRQ(ierr);
1506   ierr = HYPRE_SStructSysPFMGSetup(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x);CHKERRQ(ierr);
1507 
1508   PetscFunctionReturn(0);
1509 }
1510 
1511 
1512 /*MC
1513      PCSysPFMG - the hypre SysPFMG multigrid solver
1514 
1515    Level: advanced
1516 
1517    Options Database:
1518 + -pc_syspfmg_its <its> number of iterations of SysPFMG to use as preconditioner
1519 . -pc_syspfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
1520 . -pc_syspfmg_num_post_relax <steps> number of smoothing steps after coarse grid
1521 . -pc_syspfmg_tol <tol> tolerance of SysPFMG
1522 . -pc_syspfmg_relax_type -relaxation type for the up and down cycles, one of Weighted-Jacobi,Red/Black-Gauss-Seidel
1523 
1524    Notes:  This is for CELL-centered descretizations
1525 
1526            This must be used with the MATHYPRESSTRUCT matrix type.
1527            This is less general than in hypre, it supports only one part, and one block per process defined by a PETSc DA.
1528            Also, only cell-centered variables.
1529 
1530 .seealso:  PCMG, MATHYPRESSTRUCT
1531 M*/
1532 
1533 EXTERN_C_BEGIN
1534 #undef __FUNCT__
1535 #define __FUNCT__ "PCCreate_SysPFMG"
1536 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_SysPFMG(PC pc)
1537 {
1538   PetscErrorCode     ierr;
1539   PC_SysPFMG        *ex;
1540 
1541   PetscFunctionBegin;
1542   ierr = PetscNew(PC_SysPFMG,&ex);CHKERRQ(ierr);\
1543   pc->data = ex;
1544 
1545   ex->its            = 1;
1546   ex->tol            = 1.e-8;
1547   ex->relax_type     = 1;
1548   ex->num_pre_relax  = 1;
1549   ex->num_post_relax = 1;
1550 
1551   pc->ops->setfromoptions  = PCSetFromOptions_SysPFMG;
1552   pc->ops->view            = PCView_SysPFMG;
1553   pc->ops->destroy         = PCDestroy_SysPFMG;
1554   pc->ops->apply           = PCApply_SysPFMG;
1555   pc->ops->applyrichardson = PCApplyRichardson_SysPFMG;
1556   pc->ops->setup           = PCSetUp_SysPFMG;
1557   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr);
1558   ierr = HYPRE_SStructSysPFMGCreate(ex->hcomm,&ex->ss_solver);CHKERRQ(ierr);
1559   PetscFunctionReturn(0);
1560 }
1561 EXTERN_C_END
1562