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