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