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