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