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