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