xref: /petsc/src/ksp/pc/impls/hypre/hypre.c (revision 5f9db2b25dbeb4dbaf539bdb264ccb4c01c47bc2)
1 /*
2    Provides an interface to the LLNL package hypre
3 */
4 
5 #include <petscpkg_version.h>
6 #include <petsc/private/pcimpl.h>          /*I "petscpc.h" I*/
7 /* this include is needed ONLY to allow access to the private data inside the Mat object specific to hypre */
8 #include <petsc/private/matimpl.h>
9 #include <petsc/private/vecimpl.h>
10 #include <../src/vec/vec/impls/hypre/vhyp.h>
11 #include <../src/mat/impls/hypre/mhypre.h>
12 #include <../src/dm/impls/da/hypre/mhyp.h>
13 #include <_hypre_parcsr_ls.h>
14 #include <petscmathypre.h>
15 
16 #if defined(PETSC_HAVE_HYPRE_DEVICE)
17 #include <petsc/private/deviceimpl.h>
18 #endif
19 
20 static PetscBool cite = PETSC_FALSE;
21 static const char hypreCitation[] = "@manual{hypre-web-page,\n  title  = {{\\sl hypre}: High Performance Preconditioners},\n  organization = {Lawrence Livermore National Laboratory},\n  note  = {\\url{https://computation.llnl.gov/projects/hypre-scalable-linear-solvers-multigrid-methods}}\n}\n";
22 
23 /*
24    Private context (data structure) for the  preconditioner.
25 */
26 typedef struct {
27   HYPRE_Solver   hsolver;
28   Mat            hpmat; /* MatHYPRE */
29 
30   HYPRE_Int (*destroy)(HYPRE_Solver);
31   HYPRE_Int (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
32   HYPRE_Int (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
33 
34   MPI_Comm comm_hypre;
35   char     *hypre_type;
36 
37   /* options for Pilut and BoomerAMG*/
38   PetscInt  maxiter;
39   PetscReal tol;
40 
41   /* options for Pilut */
42   PetscInt factorrowsize;
43 
44   /* options for ParaSails */
45   PetscInt  nlevels;
46   PetscReal threshold;
47   PetscReal filter;
48   PetscReal loadbal;
49   PetscInt  logging;
50   PetscInt  ruse;
51   PetscInt  symt;
52 
53   /* options for BoomerAMG */
54   PetscBool printstatistics;
55 
56   /* options for BoomerAMG */
57   PetscInt  cycletype;
58   PetscInt  maxlevels;
59   PetscReal strongthreshold;
60   PetscReal maxrowsum;
61   PetscInt  gridsweeps[3];
62   PetscInt  coarsentype;
63   PetscInt  measuretype;
64   PetscInt  smoothtype;
65   PetscInt  smoothnumlevels;
66   PetscInt  eu_level;   /* Number of levels for ILU(k) in Euclid */
67   PetscReal eu_droptolerance; /* Drop tolerance for ILU(k) in Euclid */
68   PetscInt  eu_bj;      /* Defines use of Block Jacobi ILU in Euclid */
69   PetscInt  relaxtype[3];
70   PetscReal relaxweight;
71   PetscReal outerrelaxweight;
72   PetscInt  relaxorder;
73   PetscReal truncfactor;
74   PetscBool applyrichardson;
75   PetscInt  pmax;
76   PetscInt  interptype;
77   PetscInt  maxc;
78   PetscInt  minc;
79 #if PETSC_PKG_HYPRE_VERSION_GE(2,23,0)
80   char      *spgemm_type; // this is a global hypre parameter but is closely associated with BoomerAMG
81 #endif
82   /* GPU */
83   PetscBool keeptranspose;
84   PetscInt  rap2;
85   PetscInt  mod_rap2;
86 
87   /* AIR */
88   PetscInt  Rtype;
89   PetscReal Rstrongthreshold;
90   PetscReal Rfilterthreshold;
91   PetscInt  Adroptype;
92   PetscReal Adroptol;
93 
94   PetscInt  agg_nl;
95   PetscInt  agg_interptype;
96   PetscInt  agg_num_paths;
97   PetscBool nodal_relax;
98   PetscInt  nodal_relax_levels;
99 
100   PetscInt  nodal_coarsening;
101   PetscInt  nodal_coarsening_diag;
102   PetscInt  vec_interp_variant;
103   PetscInt  vec_interp_qmax;
104   PetscBool vec_interp_smooth;
105   PetscInt  interp_refine;
106 
107   /* NearNullSpace support */
108   VecHYPRE_IJVector *hmnull;
109   HYPRE_ParVector   *phmnull;
110   PetscInt          n_hmnull;
111   Vec               hmnull_constant;
112 
113   /* options for AS (Auxiliary Space preconditioners) */
114   PetscInt  as_print;
115   PetscInt  as_max_iter;
116   PetscReal as_tol;
117   PetscInt  as_relax_type;
118   PetscInt  as_relax_times;
119   PetscReal as_relax_weight;
120   PetscReal as_omega;
121   PetscInt  as_amg_alpha_opts[5]; /* AMG coarsen type, agg_levels, relax_type, interp_type, Pmax for vector Poisson (AMS) or Curl problem (ADS) */
122   PetscReal as_amg_alpha_theta;   /* AMG strength for vector Poisson (AMS) or Curl problem (ADS) */
123   PetscInt  as_amg_beta_opts[5];  /* AMG coarsen type, agg_levels, relax_type, interp_type, Pmax for scalar Poisson (AMS) or vector Poisson (ADS) */
124   PetscReal as_amg_beta_theta;    /* AMG strength for scalar Poisson (AMS) or vector Poisson (ADS)  */
125   PetscInt  ams_cycle_type;
126   PetscInt  ads_cycle_type;
127 
128   /* additional data */
129   Mat G;             /* MatHYPRE */
130   Mat C;             /* MatHYPRE */
131   Mat alpha_Poisson; /* MatHYPRE */
132   Mat beta_Poisson;  /* MatHYPRE */
133 
134   /* extra information for AMS */
135   PetscInt          dim; /* geometrical dimension */
136   VecHYPRE_IJVector coords[3];
137   VecHYPRE_IJVector constants[3];
138   VecHYPRE_IJVector interior;
139   Mat               RT_PiFull, RT_Pi[3];
140   Mat               ND_PiFull, ND_Pi[3];
141   PetscBool         ams_beta_is_zero;
142   PetscBool         ams_beta_is_zero_part;
143   PetscInt          ams_proj_freq;
144 } PC_HYPRE;
145 
146 PetscErrorCode PCHYPREGetSolver(PC pc,HYPRE_Solver *hsolver)
147 {
148   PC_HYPRE *jac = (PC_HYPRE*)pc->data;
149 
150   PetscFunctionBegin;
151   *hsolver = jac->hsolver;
152   PetscFunctionReturn(0);
153 }
154 
155 /*
156   Matrices with AIJ format are created IN PLACE with using (I,J,data) from BoomerAMG. Since the data format in hypre_ParCSRMatrix
157   is different from that used in PETSc, the original hypre_ParCSRMatrix can not be used any more after call this routine.
158   It is used in PCHMG. Other users should avoid using this function.
159 */
160 static PetscErrorCode PCGetCoarseOperators_BoomerAMG(PC pc,PetscInt *nlevels,Mat *operators[])
161 {
162   PC_HYPRE             *jac  = (PC_HYPRE*)pc->data;
163   PetscBool            same = PETSC_FALSE;
164   PetscInt             num_levels,l;
165   Mat                  *mattmp;
166   hypre_ParCSRMatrix   **A_array;
167 
168   PetscFunctionBegin;
169   PetscCall(PetscStrcmp(jac->hypre_type,"boomeramg",&same));
170   PetscCheck(same,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_NOTSAMETYPE,"Hypre type is not BoomerAMG ");
171   num_levels = hypre_ParAMGDataNumLevels((hypre_ParAMGData*) (jac->hsolver));
172   PetscCall(PetscMalloc1(num_levels,&mattmp));
173   A_array    = hypre_ParAMGDataAArray((hypre_ParAMGData*) (jac->hsolver));
174   for (l=1; l<num_levels; l++) {
175     PetscCall(MatCreateFromParCSR(A_array[l],MATAIJ,PETSC_OWN_POINTER, &(mattmp[num_levels-1-l])));
176     /* We want to own the data, and HYPRE can not touch this matrix any more */
177     A_array[l] = NULL;
178   }
179   *nlevels = num_levels;
180   *operators = mattmp;
181   PetscFunctionReturn(0);
182 }
183 
184 /*
185   Matrices with AIJ format are created IN PLACE with using (I,J,data) from BoomerAMG. Since the data format in hypre_ParCSRMatrix
186   is different from that used in PETSc, the original hypre_ParCSRMatrix can not be used any more after call this routine.
187   It is used in PCHMG. Other users should avoid using this function.
188 */
189 static PetscErrorCode PCGetInterpolations_BoomerAMG(PC pc,PetscInt *nlevels,Mat *interpolations[])
190 {
191   PC_HYPRE             *jac  = (PC_HYPRE*)pc->data;
192   PetscBool            same = PETSC_FALSE;
193   PetscInt             num_levels,l;
194   Mat                  *mattmp;
195   hypre_ParCSRMatrix   **P_array;
196 
197   PetscFunctionBegin;
198   PetscCall(PetscStrcmp(jac->hypre_type,"boomeramg",&same));
199   PetscCheck(same,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_NOTSAMETYPE,"Hypre type is not BoomerAMG ");
200   num_levels = hypre_ParAMGDataNumLevels((hypre_ParAMGData*) (jac->hsolver));
201   PetscCall(PetscMalloc1(num_levels,&mattmp));
202   P_array  = hypre_ParAMGDataPArray((hypre_ParAMGData*) (jac->hsolver));
203   for (l=1; l<num_levels; l++) {
204     PetscCall(MatCreateFromParCSR(P_array[num_levels-1-l],MATAIJ,PETSC_OWN_POINTER, &(mattmp[l-1])));
205     /* We want to own the data, and HYPRE can not touch this matrix any more */
206     P_array[num_levels-1-l] = NULL;
207   }
208   *nlevels = num_levels;
209   *interpolations = mattmp;
210   PetscFunctionReturn(0);
211 }
212 
213 /* Resets (frees) Hypre's representation of the near null space */
214 static PetscErrorCode PCHYPREResetNearNullSpace_Private(PC pc)
215 {
216   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
217   PetscInt       i;
218 
219   PetscFunctionBegin;
220   for (i=0; i<jac->n_hmnull; i++) {
221     PetscCall(VecHYPRE_IJVectorDestroy(&jac->hmnull[i]));
222   }
223   PetscCall(PetscFree(jac->hmnull));
224   PetscCall(PetscFree(jac->phmnull));
225   PetscCall(VecDestroy(&jac->hmnull_constant));
226   jac->n_hmnull = 0;
227   PetscFunctionReturn(0);
228 }
229 
230 static PetscErrorCode PCSetUp_HYPRE(PC pc)
231 {
232   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
233   Mat_HYPRE          *hjac;
234   HYPRE_ParCSRMatrix hmat;
235   HYPRE_ParVector    bv,xv;
236   PetscBool          ishypre;
237 
238   PetscFunctionBegin;
239   if (!jac->hypre_type) {
240     PetscCall(PCHYPRESetType(pc,"boomeramg"));
241   }
242 
243   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRE,&ishypre));
244   if (!ishypre) {
245     PetscCall(MatDestroy(&jac->hpmat));
246     PetscCall(MatConvert(pc->pmat,MATHYPRE,MAT_INITIAL_MATRIX,&jac->hpmat));
247     PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->hpmat));
248   } else {
249     PetscCall(PetscObjectReference((PetscObject)pc->pmat));
250     PetscCall(MatDestroy(&jac->hpmat));
251     jac->hpmat = pc->pmat;
252   }
253   /* allow debug */
254   PetscCall(MatViewFromOptions(jac->hpmat,NULL,"-pc_hypre_mat_view"));
255   hjac = (Mat_HYPRE*)(jac->hpmat->data);
256 
257   /* special case for BoomerAMG */
258   if (jac->setup == HYPRE_BoomerAMGSetup) {
259     MatNullSpace mnull;
260     PetscBool    has_const;
261     PetscInt     bs,nvec,i;
262     const Vec    *vecs;
263 
264     PetscCall(MatGetBlockSize(pc->pmat,&bs));
265     if (bs > 1) PetscCallExternal(HYPRE_BoomerAMGSetNumFunctions,jac->hsolver,bs);
266     PetscCall(MatGetNearNullSpace(pc->mat, &mnull));
267     if (mnull) {
268       PetscCall(PCHYPREResetNearNullSpace_Private(pc));
269       PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
270       PetscCall(PetscMalloc1(nvec+1,&jac->hmnull));
271       PetscCall(PetscMalloc1(nvec+1,&jac->phmnull));
272       for (i=0; i<nvec; i++) {
273         PetscCall(VecHYPRE_IJVectorCreate(vecs[i]->map,&jac->hmnull[i]));
274         PetscCall(VecHYPRE_IJVectorCopy(vecs[i],jac->hmnull[i]));
275         PetscCallExternal(HYPRE_IJVectorGetObject,jac->hmnull[i]->ij,(void**)&jac->phmnull[i]);
276       }
277       if (has_const) {
278         PetscCall(MatCreateVecs(pc->pmat,&jac->hmnull_constant,NULL));
279         PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)jac->hmnull_constant));
280         PetscCall(VecSet(jac->hmnull_constant,1));
281         PetscCall(VecNormalize(jac->hmnull_constant,NULL));
282         PetscCall(VecHYPRE_IJVectorCreate(jac->hmnull_constant->map,&jac->hmnull[nvec]));
283         PetscCall(VecHYPRE_IJVectorCopy(jac->hmnull_constant,jac->hmnull[nvec]));
284         PetscCallExternal(HYPRE_IJVectorGetObject,jac->hmnull[nvec]->ij,(void**)&jac->phmnull[nvec]);
285         nvec++;
286       }
287       PetscCallExternal(HYPRE_BoomerAMGSetInterpVectors,jac->hsolver,nvec,jac->phmnull);
288       jac->n_hmnull = nvec;
289     }
290   }
291 
292   /* special case for AMS */
293   if (jac->setup == HYPRE_AMSSetup) {
294     Mat_HYPRE          *hm;
295     HYPRE_ParCSRMatrix parcsr;
296     if (!jac->coords[0] && !jac->constants[0] && !(jac->ND_PiFull || (jac->ND_Pi[0] && jac->ND_Pi[1]))) {
297       SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE AMS preconditioner needs either the coordinate vectors via PCSetCoordinates() or the edge constant vectors via PCHYPRESetEdgeConstantVectors() or the interpolation matrix via PCHYPRESetInterpolations");
298     }
299     if (jac->dim) {
300       PetscCallExternal(HYPRE_AMSSetDimension,jac->hsolver,jac->dim);
301     }
302     if (jac->constants[0]) {
303       HYPRE_ParVector ozz,zoz,zzo = NULL;
304       PetscCallExternal(HYPRE_IJVectorGetObject,jac->constants[0]->ij,(void**)(&ozz));
305       PetscCallExternal(HYPRE_IJVectorGetObject,jac->constants[1]->ij,(void**)(&zoz));
306       if (jac->constants[2]) {
307         PetscCallExternal(HYPRE_IJVectorGetObject,jac->constants[2]->ij,(void**)(&zzo));
308       }
309       PetscCallExternal(HYPRE_AMSSetEdgeConstantVectors,jac->hsolver,ozz,zoz,zzo);
310     }
311     if (jac->coords[0]) {
312       HYPRE_ParVector coords[3];
313       coords[0] = NULL;
314       coords[1] = NULL;
315       coords[2] = NULL;
316       if (jac->coords[0]) PetscCallExternal(HYPRE_IJVectorGetObject,jac->coords[0]->ij,(void**)(&coords[0]));
317       if (jac->coords[1]) PetscCallExternal(HYPRE_IJVectorGetObject,jac->coords[1]->ij,(void**)(&coords[1]));
318       if (jac->coords[2]) PetscCallExternal(HYPRE_IJVectorGetObject,jac->coords[2]->ij,(void**)(&coords[2]));
319       PetscCallExternal(HYPRE_AMSSetCoordinateVectors,jac->hsolver,coords[0],coords[1],coords[2]);
320     }
321     PetscCheck(jac->G,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE AMS preconditioner needs the discrete gradient operator via PCHYPRESetDiscreteGradient");
322     hm = (Mat_HYPRE*)(jac->G->data);
323     PetscCallExternal(HYPRE_IJMatrixGetObject,hm->ij,(void**)(&parcsr));
324     PetscCallExternal(HYPRE_AMSSetDiscreteGradient,jac->hsolver,parcsr);
325     if (jac->alpha_Poisson) {
326       hm = (Mat_HYPRE*)(jac->alpha_Poisson->data);
327       PetscCallExternal(HYPRE_IJMatrixGetObject,hm->ij,(void**)(&parcsr));
328       PetscCallExternal(HYPRE_AMSSetAlphaPoissonMatrix,jac->hsolver,parcsr);
329     }
330     if (jac->ams_beta_is_zero) {
331       PetscCallExternal(HYPRE_AMSSetBetaPoissonMatrix,jac->hsolver,NULL);
332     } else if (jac->beta_Poisson) {
333       hm = (Mat_HYPRE*)(jac->beta_Poisson->data);
334       PetscCallExternal(HYPRE_IJMatrixGetObject,hm->ij,(void**)(&parcsr));
335       PetscCallExternal(HYPRE_AMSSetBetaPoissonMatrix,jac->hsolver,parcsr);
336     } else if (jac->ams_beta_is_zero_part)  {
337       if (jac->interior) {
338         HYPRE_ParVector interior = NULL;
339         PetscCallExternal(HYPRE_IJVectorGetObject,jac->interior->ij,(void**)(&interior));
340         PetscCallExternal(HYPRE_AMSSetInteriorNodes,jac->hsolver,interior);
341       } else {
342         jac->ams_beta_is_zero_part = PETSC_FALSE;
343       }
344     }
345     if (jac->ND_PiFull || (jac->ND_Pi[0] && jac->ND_Pi[1])) {
346       PetscInt           i;
347       HYPRE_ParCSRMatrix nd_parcsrfull, nd_parcsr[3];
348       if (jac->ND_PiFull) {
349         hm = (Mat_HYPRE*)(jac->ND_PiFull->data);
350         PetscCallExternal(HYPRE_IJMatrixGetObject,hm->ij,(void**)(&nd_parcsrfull));
351       } else {
352         nd_parcsrfull = NULL;
353       }
354       for (i=0;i<3;++i) {
355         if (jac->ND_Pi[i]) {
356           hm = (Mat_HYPRE*)(jac->ND_Pi[i]->data);
357           PetscCallExternal(HYPRE_IJMatrixGetObject,hm->ij,(void**)(&nd_parcsr[i]));
358         } else {
359           nd_parcsr[i] = NULL;
360         }
361       }
362       PetscCallExternal(HYPRE_AMSSetInterpolations,jac->hsolver,nd_parcsrfull,nd_parcsr[0],nd_parcsr[1],nd_parcsr[2]);
363     }
364   }
365   /* special case for ADS */
366   if (jac->setup == HYPRE_ADSSetup) {
367     Mat_HYPRE          *hm;
368     HYPRE_ParCSRMatrix parcsr;
369     if (!jac->coords[0] && !((jac->RT_PiFull || (jac->RT_Pi[0] && jac->RT_Pi[1])) && (jac->ND_PiFull || (jac->ND_Pi[0] && jac->ND_Pi[1])))) {
370       SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE ADS preconditioner needs either the coordinate vectors via PCSetCoordinates() or the interpolation matrices via PCHYPRESetInterpolations");
371     }
372     else PetscCheck(jac->coords[1] && jac->coords[2],PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE ADS preconditioner has been designed for three dimensional problems! For two dimensional problems, use HYPRE AMS instead");
373     PetscCheck(jac->G,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE ADS preconditioner needs the discrete gradient operator via PCHYPRESetDiscreteGradient");
374     PetscCheck(jac->C,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE ADS preconditioner needs the discrete curl operator via PCHYPRESetDiscreteGradient");
375     if (jac->coords[0]) {
376       HYPRE_ParVector coords[3];
377       coords[0] = NULL;
378       coords[1] = NULL;
379       coords[2] = NULL;
380       if (jac->coords[0]) PetscCallExternal(HYPRE_IJVectorGetObject,jac->coords[0]->ij,(void**)(&coords[0]));
381       if (jac->coords[1]) PetscCallExternal(HYPRE_IJVectorGetObject,jac->coords[1]->ij,(void**)(&coords[1]));
382       if (jac->coords[2]) PetscCallExternal(HYPRE_IJVectorGetObject,jac->coords[2]->ij,(void**)(&coords[2]));
383       PetscCallExternal(HYPRE_ADSSetCoordinateVectors,jac->hsolver,coords[0],coords[1],coords[2]);
384     }
385     hm = (Mat_HYPRE*)(jac->G->data);
386     PetscCallExternal(HYPRE_IJMatrixGetObject,hm->ij,(void**)(&parcsr));
387     PetscCallExternal(HYPRE_ADSSetDiscreteGradient,jac->hsolver,parcsr);
388     hm = (Mat_HYPRE*)(jac->C->data);
389     PetscCallExternal(HYPRE_IJMatrixGetObject,hm->ij,(void**)(&parcsr));
390     PetscCallExternal(HYPRE_ADSSetDiscreteCurl,jac->hsolver,parcsr);
391     if ((jac->RT_PiFull || (jac->RT_Pi[0] && jac->RT_Pi[1])) && (jac->ND_PiFull || (jac->ND_Pi[0] && jac->ND_Pi[1]))) {
392       PetscInt           i;
393       HYPRE_ParCSRMatrix rt_parcsrfull, rt_parcsr[3];
394       HYPRE_ParCSRMatrix nd_parcsrfull, nd_parcsr[3];
395       if (jac->RT_PiFull) {
396         hm = (Mat_HYPRE*)(jac->RT_PiFull->data);
397         PetscCallExternal(HYPRE_IJMatrixGetObject,hm->ij,(void**)(&rt_parcsrfull));
398       } else {
399         rt_parcsrfull = NULL;
400       }
401       for (i=0;i<3;++i) {
402         if (jac->RT_Pi[i]) {
403           hm = (Mat_HYPRE*)(jac->RT_Pi[i]->data);
404           PetscCallExternal(HYPRE_IJMatrixGetObject,hm->ij,(void**)(&rt_parcsr[i]));
405         } else {
406           rt_parcsr[i] = NULL;
407         }
408       }
409       if (jac->ND_PiFull) {
410         hm = (Mat_HYPRE*)(jac->ND_PiFull->data);
411         PetscCallExternal(HYPRE_IJMatrixGetObject,hm->ij,(void**)(&nd_parcsrfull));
412       } else {
413         nd_parcsrfull = NULL;
414       }
415       for (i=0;i<3;++i) {
416         if (jac->ND_Pi[i]) {
417           hm = (Mat_HYPRE*)(jac->ND_Pi[i]->data);
418           PetscCallExternal(HYPRE_IJMatrixGetObject,hm->ij,(void**)(&nd_parcsr[i]));
419         } else {
420           nd_parcsr[i] = NULL;
421         }
422       }
423       PetscCallExternal(HYPRE_ADSSetInterpolations,jac->hsolver,rt_parcsrfull,rt_parcsr[0],rt_parcsr[1],rt_parcsr[2],nd_parcsrfull,nd_parcsr[0],nd_parcsr[1],nd_parcsr[2]);
424     }
425   }
426   PetscCallExternal(HYPRE_IJMatrixGetObject,hjac->ij,(void**)&hmat);
427   PetscCallExternal(HYPRE_IJVectorGetObject,hjac->b->ij,(void**)&bv);
428   PetscCallExternal(HYPRE_IJVectorGetObject,hjac->x->ij,(void**)&xv);
429   PetscCallExternal(jac->setup,jac->hsolver,hmat,bv,xv);
430   PetscFunctionReturn(0);
431 }
432 
433 static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
434 {
435   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
436   Mat_HYPRE          *hjac = (Mat_HYPRE*)(jac->hpmat->data);
437   HYPRE_ParCSRMatrix hmat;
438   HYPRE_ParVector    jbv,jxv;
439 
440   PetscFunctionBegin;
441   PetscCall(PetscCitationsRegister(hypreCitation,&cite));
442   if (!jac->applyrichardson) PetscCall(VecSet(x,0.0));
443   PetscCall(VecHYPRE_IJVectorPushVecRead(hjac->b,b));
444   if (jac->applyrichardson) PetscCall(VecHYPRE_IJVectorPushVec(hjac->x,x));
445   else PetscCall(VecHYPRE_IJVectorPushVecWrite(hjac->x,x));
446   PetscCallExternal(HYPRE_IJMatrixGetObject,hjac->ij,(void**)&hmat);
447   PetscCallExternal(HYPRE_IJVectorGetObject,hjac->b->ij,(void**)&jbv);
448   PetscCallExternal(HYPRE_IJVectorGetObject,hjac->x->ij,(void**)&jxv);
449   PetscStackCallExternalVoid("Hypre solve",do {
450       HYPRE_Int hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv);
451       if (hierr) {
452         PetscCheck(hierr == HYPRE_ERROR_CONV,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",(int)hierr);
453         hypre__global_error = 0;
454       }
455     } while (0));
456 
457   if (jac->setup == HYPRE_AMSSetup && jac->ams_beta_is_zero_part) {
458     PetscCallExternal(HYPRE_AMSProjectOutGradients,jac->hsolver,jxv);
459   }
460   PetscCall(VecHYPRE_IJVectorPopVec(hjac->x));
461   PetscCall(VecHYPRE_IJVectorPopVec(hjac->b));
462   PetscFunctionReturn(0);
463 }
464 
465 static PetscErrorCode PCReset_HYPRE(PC pc)
466 {
467   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
468 
469   PetscFunctionBegin;
470   PetscCall(MatDestroy(&jac->hpmat));
471   PetscCall(MatDestroy(&jac->G));
472   PetscCall(MatDestroy(&jac->C));
473   PetscCall(MatDestroy(&jac->alpha_Poisson));
474   PetscCall(MatDestroy(&jac->beta_Poisson));
475   PetscCall(MatDestroy(&jac->RT_PiFull));
476   PetscCall(MatDestroy(&jac->RT_Pi[0]));
477   PetscCall(MatDestroy(&jac->RT_Pi[1]));
478   PetscCall(MatDestroy(&jac->RT_Pi[2]));
479   PetscCall(MatDestroy(&jac->ND_PiFull));
480   PetscCall(MatDestroy(&jac->ND_Pi[0]));
481   PetscCall(MatDestroy(&jac->ND_Pi[1]));
482   PetscCall(MatDestroy(&jac->ND_Pi[2]));
483   PetscCall(VecHYPRE_IJVectorDestroy(&jac->coords[0]));
484   PetscCall(VecHYPRE_IJVectorDestroy(&jac->coords[1]));
485   PetscCall(VecHYPRE_IJVectorDestroy(&jac->coords[2]));
486   PetscCall(VecHYPRE_IJVectorDestroy(&jac->constants[0]));
487   PetscCall(VecHYPRE_IJVectorDestroy(&jac->constants[1]));
488   PetscCall(VecHYPRE_IJVectorDestroy(&jac->constants[2]));
489   PetscCall(VecHYPRE_IJVectorDestroy(&jac->interior));
490   PetscCall(PCHYPREResetNearNullSpace_Private(pc));
491   jac->ams_beta_is_zero = PETSC_FALSE;
492   jac->ams_beta_is_zero_part = PETSC_FALSE;
493   jac->dim = 0;
494   PetscFunctionReturn(0);
495 }
496 
497 static PetscErrorCode PCDestroy_HYPRE(PC pc)
498 {
499   PC_HYPRE                 *jac = (PC_HYPRE*)pc->data;
500 
501   PetscFunctionBegin;
502   PetscCall(PCReset_HYPRE(pc));
503   if (jac->destroy) PetscCallExternal(jac->destroy,jac->hsolver);
504   PetscCall(PetscFree(jac->hypre_type));
505 #if PETSC_PKG_HYPRE_VERSION_GE(2,23,0)
506   PetscCall(PetscFree(jac->spgemm_type));
507 #endif
508   if (jac->comm_hypre != MPI_COMM_NULL) PetscCall(PetscCommRestoreComm(PetscObjectComm((PetscObject)pc),&jac->comm_hypre));
509   PetscCall(PetscFree(pc->data));
510 
511   PetscCall(PetscObjectChangeTypeName((PetscObject)pc,0));
512   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetType_C",NULL));
513   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPREGetType_C",NULL));
514   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteGradient_C",NULL));
515   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteCurl_C",NULL));
516   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetInterpolations_C",NULL));
517   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetConstantEdgeVectors_C",NULL));
518   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetPoissonMatrix_C",NULL));
519   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetEdgeConstantVectors_C",NULL));
520   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPREAMSSetInteriorNodes_C",NULL));
521   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL));
522   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL));
523   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGalerkinSetMatProductAlgorithm_C",NULL));
524   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGalerkinGetMatProductAlgorithm_C",NULL));
525   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",NULL));
526   PetscFunctionReturn(0);
527 }
528 
529 /* --------------------------------------------------------------------------------------------*/
530 static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PetscOptionItems *PetscOptionsObject,PC pc)
531 {
532   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
533   PetscBool      flag;
534 
535   PetscFunctionBegin;
536   PetscOptionsHeadBegin(PetscOptionsObject,"HYPRE Pilut Options");
537   PetscCall(PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag));
538   if (flag) PetscCallExternal(HYPRE_ParCSRPilutSetMaxIter,jac->hsolver,jac->maxiter);
539   PetscCall(PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag));
540   if (flag) PetscCallExternal(HYPRE_ParCSRPilutSetDropTolerance,jac->hsolver,jac->tol);
541   PetscCall(PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag));
542   if (flag) PetscCallExternal(HYPRE_ParCSRPilutSetFactorRowSize,jac->hsolver,jac->factorrowsize);
543   PetscOptionsHeadEnd();
544   PetscFunctionReturn(0);
545 }
546 
547 static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
548 {
549   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
550   PetscBool      iascii;
551 
552   PetscFunctionBegin;
553   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
554   if (iascii) {
555     PetscCall(PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n"));
556     if (jac->maxiter != PETSC_DEFAULT) {
557       PetscCall(PetscViewerASCIIPrintf(viewer,"    maximum number of iterations %" PetscInt_FMT "\n",jac->maxiter));
558     } else {
559       PetscCall(PetscViewerASCIIPrintf(viewer,"    default maximum number of iterations \n"));
560     }
561     if (jac->tol != PETSC_DEFAULT) {
562       PetscCall(PetscViewerASCIIPrintf(viewer,"    drop tolerance %g\n",(double)jac->tol));
563     } else {
564       PetscCall(PetscViewerASCIIPrintf(viewer,"    default drop tolerance \n"));
565     }
566     if (jac->factorrowsize != PETSC_DEFAULT) {
567       PetscCall(PetscViewerASCIIPrintf(viewer,"    factor row size %" PetscInt_FMT "\n",jac->factorrowsize));
568     } else {
569       PetscCall(PetscViewerASCIIPrintf(viewer,"    default factor row size \n"));
570     }
571   }
572   PetscFunctionReturn(0);
573 }
574 
575 /* --------------------------------------------------------------------------------------------*/
576 static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PetscOptionItems *PetscOptionsObject,PC pc)
577 {
578   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
579   PetscBool      flag,eu_bj = jac->eu_bj ? PETSC_TRUE : PETSC_FALSE;
580 
581   PetscFunctionBegin;
582   PetscOptionsHeadBegin(PetscOptionsObject,"HYPRE Euclid Options");
583   PetscCall(PetscOptionsInt("-pc_hypre_euclid_level","Factorization levels","None",jac->eu_level,&jac->eu_level,&flag));
584   if (flag) PetscCallExternal(HYPRE_EuclidSetLevel,jac->hsolver,jac->eu_level);
585 
586   PetscCall(PetscOptionsReal("-pc_hypre_euclid_droptolerance","Drop tolerance for ILU(k) in Euclid","None",jac->eu_droptolerance,&jac->eu_droptolerance,&flag));
587   if (flag) {
588     PetscMPIInt size;
589 
590     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size));
591     PetscCheck(size == 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"hypre's Euclid does not support a parallel drop tolerance");
592     PetscCallExternal(HYPRE_EuclidSetILUT,jac->hsolver,jac->eu_droptolerance);
593   }
594 
595   PetscCall(PetscOptionsBool("-pc_hypre_euclid_bj", "Use Block Jacobi for ILU in Euclid", "None", eu_bj,&eu_bj,&flag));
596   if (flag) {
597     jac->eu_bj = eu_bj ? 1 : 0;
598     PetscCallExternal(HYPRE_EuclidSetBJ,jac->hsolver,jac->eu_bj);
599   }
600   PetscOptionsHeadEnd();
601   PetscFunctionReturn(0);
602 }
603 
604 static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
605 {
606   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
607   PetscBool      iascii;
608 
609   PetscFunctionBegin;
610   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
611   if (iascii) {
612     PetscCall(PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n"));
613     if (jac->eu_level != PETSC_DEFAULT) {
614       PetscCall(PetscViewerASCIIPrintf(viewer,"    factorization levels %" PetscInt_FMT "\n",jac->eu_level));
615     } else {
616       PetscCall(PetscViewerASCIIPrintf(viewer,"    default factorization levels \n"));
617     }
618     PetscCall(PetscViewerASCIIPrintf(viewer,"    drop tolerance %g\n",(double)jac->eu_droptolerance));
619     PetscCall(PetscViewerASCIIPrintf(viewer,"    use Block-Jacobi? %" PetscInt_FMT "\n",jac->eu_bj));
620   }
621   PetscFunctionReturn(0);
622 }
623 
624 /* --------------------------------------------------------------------------------------------*/
625 
626 static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
627 {
628   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
629   Mat_HYPRE          *hjac = (Mat_HYPRE*)(jac->hpmat->data);
630   HYPRE_ParCSRMatrix hmat;
631   HYPRE_ParVector    jbv,jxv;
632 
633   PetscFunctionBegin;
634   PetscCall(PetscCitationsRegister(hypreCitation,&cite));
635   PetscCall(VecSet(x,0.0));
636   PetscCall(VecHYPRE_IJVectorPushVecRead(hjac->x,b));
637   PetscCall(VecHYPRE_IJVectorPushVecWrite(hjac->b,x));
638 
639   PetscCallExternal(HYPRE_IJMatrixGetObject,hjac->ij,(void**)&hmat);
640   PetscCallExternal(HYPRE_IJVectorGetObject,hjac->b->ij,(void**)&jbv);
641   PetscCallExternal(HYPRE_IJVectorGetObject,hjac->x->ij,(void**)&jxv);
642 
643   PetscStackCallExternalVoid("Hypre Transpose solve",do {
644       HYPRE_Int hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv);
645       if (hierr) {
646         /* error code of 1 in BoomerAMG merely means convergence not achieved */
647         PetscCheck(hierr == 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",(int)hierr);
648         hypre__global_error = 0;
649       }
650     } while (0));
651 
652   PetscCall(VecHYPRE_IJVectorPopVec(hjac->x));
653   PetscCall(VecHYPRE_IJVectorPopVec(hjac->b));
654   PetscFunctionReturn(0);
655 }
656 
657 static PetscErrorCode PCMGGalerkinSetMatProductAlgorithm_HYPRE_BoomerAMG(PC pc,const char name[])
658 {
659   PC_HYPRE *jac  = (PC_HYPRE*)pc->data;
660   PetscBool      flag;
661 
662 #if PETSC_PKG_HYPRE_VERSION_GE(2,23,0)
663   PetscFunctionBegin;
664   if (jac->spgemm_type) {
665     PetscCall(PetscStrcmp(jac->spgemm_type,name,&flag));
666     PetscCheck(flag,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Cannot reset the HYPRE SpGEMM (really we can)");
667     PetscFunctionReturn(0);
668   } else {
669     PetscCall(PetscStrallocpy(name, &jac->spgemm_type));
670   }
671   PetscCall(PetscStrcmp("cusparse",jac->spgemm_type,&flag));
672   if (flag) {
673     PetscCallExternal(HYPRE_SetSpGemmUseCusparse,1);
674     PetscFunctionReturn(0);
675   }
676   PetscCall(PetscStrcmp("hypre",jac->spgemm_type,&flag));
677   if (flag) {
678     PetscCallExternal(HYPRE_SetSpGemmUseCusparse,0);
679     PetscFunctionReturn(0);
680   }
681   jac->spgemm_type = NULL;
682   SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE SpGEM type %s; Choices are cusparse, hypre",name);
683 #endif
684 }
685 
686 static PetscErrorCode PCMGGalerkinGetMatProductAlgorithm_HYPRE_BoomerAMG(PC pc, const char *spgemm[])
687 {
688   PC_HYPRE *jac  = (PC_HYPRE*)pc->data;
689 
690   PetscFunctionBegin;
691   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
692 #if PETSC_PKG_HYPRE_VERSION_GE(2,23,0)
693   *spgemm = jac->spgemm_type;
694 #endif
695   PetscFunctionReturn(0);
696 }
697 
698 static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
699 static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"};
700 static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
701 /* The following corresponds to HYPRE_BoomerAMGSetRelaxType which has many missing numbers in the enum */
702 static const char *HYPREBoomerAMGSmoothType[]  = {"Schwarz-smoothers","Pilut","ParaSails","Euclid"};
703 static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","seqboundary-Gauss-Seidel","SOR/Jacobi","backward-SOR/Jacobi",
704                                                   "" /* [5] hybrid chaotic Gauss-Seidel (works only with OpenMP) */,"symmetric-SOR/Jacobi",
705                                                   "" /* 7 */,"l1scaled-SOR/Jacobi","Gaussian-elimination",
706                                                   "" /* 10 */, "" /* 11 */, "" /* 12 */, "l1-Gauss-Seidel" /* nonsymmetric */, "backward-l1-Gauss-Seidel" /* nonsymmetric */,
707                                                   "CG" /* non-stationary */,"Chebyshev","FCF-Jacobi","l1scaled-Jacobi"};
708 static const char *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i",
709                                                   "ext+i-cc", "standard", "standard-wts", "block", "block-wtd", "FF", "FF1",
710                                                   "ext", "ad-wts", "ext-mm", "ext+i-mm", "ext+e-mm"};
711 static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PetscOptionItems *PetscOptionsObject,PC pc)
712 {
713   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
714   PetscInt       bs,n,indx,level;
715   PetscBool      flg, tmp_truth;
716   double         tmpdbl, twodbl[2];
717   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};
718   const char     *PCHYPRESpgemmTypes[] = {"cusparse","hypre"};
719 
720   PetscFunctionBegin;
721   PetscOptionsHeadBegin(PetscOptionsObject,"HYPRE BoomerAMG Options");
722   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType+1,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg));
723   if (flg) {
724     jac->cycletype = indx+1;
725     PetscCallExternal(HYPRE_BoomerAMGSetCycleType,jac->hsolver,jac->cycletype);
726   }
727   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg));
728   if (flg) {
729     PetscCheck(jac->maxlevels >= 2,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %" PetscInt_FMT " must be at least two",jac->maxlevels);
730     PetscCallExternal(HYPRE_BoomerAMGSetMaxLevels,jac->hsolver,jac->maxlevels);
731   }
732   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg));
733   if (flg) {
734     PetscCheck(jac->maxiter >= 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %" PetscInt_FMT " must be at least one",jac->maxiter);
735     PetscCallExternal(HYPRE_BoomerAMGSetMaxIter,jac->hsolver,jac->maxiter);
736   }
737   PetscCall(PetscOptionsReal("-pc_hypre_boomeramg_tol","Convergence tolerance PER hypre call (0.0 = use a fixed number of iterations)","None",jac->tol,&jac->tol,&flg));
738   if (flg) {
739     PetscCheck(jac->tol >= 0.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %g must be greater than or equal to zero",(double)jac->tol);
740     PetscCallExternal(HYPRE_BoomerAMGSetTol,jac->hsolver,jac->tol);
741   }
742   bs = 1;
743   if (pc->pmat) {
744     PetscCall(MatGetBlockSize(pc->pmat,&bs));
745   }
746   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_numfunctions","Number of functions","HYPRE_BoomerAMGSetNumFunctions",bs,&bs,&flg));
747   if (flg) {
748     PetscCallExternal(HYPRE_BoomerAMGSetNumFunctions,jac->hsolver,bs);
749   }
750 
751   PetscCall(PetscOptionsReal("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg));
752   if (flg) {
753     PetscCheck(jac->truncfactor >= 0.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %g must be great than or equal zero",(double)jac->truncfactor);
754     PetscCallExternal(HYPRE_BoomerAMGSetTruncFactor,jac->hsolver,jac->truncfactor);
755   }
756 
757   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_P_max","Max elements per row for interpolation operator (0=unlimited)","None",jac->pmax,&jac->pmax,&flg));
758   if (flg) {
759     PetscCheck(jac->pmax >= 0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"P_max %" PetscInt_FMT " must be greater than or equal to zero",jac->pmax);
760     PetscCallExternal(HYPRE_BoomerAMGSetPMaxElmts,jac->hsolver,jac->pmax);
761   }
762 
763   PetscCall(PetscOptionsRangeInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg,0,jac->maxlevels));
764   if (flg) PetscCallExternal(HYPRE_BoomerAMGSetAggNumLevels,jac->hsolver,jac->agg_nl);
765 
766   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_agg_num_paths","Number of paths for aggressive coarsening","None",jac->agg_num_paths,&jac->agg_num_paths,&flg));
767   if (flg) {
768     PetscCheck(jac->agg_num_paths >= 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of paths %" PetscInt_FMT " must be greater than or equal to 1",jac->agg_num_paths);
769     PetscCallExternal(HYPRE_BoomerAMGSetNumPaths,jac->hsolver,jac->agg_num_paths);
770   }
771 
772   PetscCall(PetscOptionsReal("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg));
773   if (flg) {
774     PetscCheck(jac->strongthreshold >= 0.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %g must be great than or equal zero",(double)jac->strongthreshold);
775     PetscCallExternal(HYPRE_BoomerAMGSetStrongThreshold,jac->hsolver,jac->strongthreshold);
776   }
777   PetscCall(PetscOptionsReal("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg));
778   if (flg) {
779     PetscCheck(jac->maxrowsum >= 0.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %g must be greater than zero",(double)jac->maxrowsum);
780     PetscCheck(jac->maxrowsum <= 1.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %g must be less than or equal one",(double)jac->maxrowsum);
781     PetscCallExternal(HYPRE_BoomerAMGSetMaxRowSum,jac->hsolver,jac->maxrowsum);
782   }
783 
784   /* Grid sweeps */
785   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_all","Number of sweeps for the up and down grid levels","None",jac->gridsweeps[0],&indx,&flg));
786   if (flg) {
787     PetscCallExternal(HYPRE_BoomerAMGSetNumSweeps,jac->hsolver,indx);
788     /* modify the jac structure so we can view the updated options with PC_View */
789     jac->gridsweeps[0] = indx;
790     jac->gridsweeps[1] = indx;
791     /*defaults coarse to 1 */
792     jac->gridsweeps[2] = 1;
793   }
794   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_nodal_coarsen","Use a nodal based coarsening 1-6","HYPRE_BoomerAMGSetNodal",jac->nodal_coarsening,&jac->nodal_coarsening,&flg));
795   if (flg) {
796     PetscCallExternal(HYPRE_BoomerAMGSetNodal,jac->hsolver,jac->nodal_coarsening);
797   }
798   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_nodal_coarsen_diag","Diagonal in strength matrix for nodal based coarsening 0-2","HYPRE_BoomerAMGSetNodalDiag",jac->nodal_coarsening_diag,&jac->nodal_coarsening_diag,&flg));
799   if (flg) {
800     PetscCallExternal(HYPRE_BoomerAMGSetNodalDiag,jac->hsolver,jac->nodal_coarsening_diag);
801   }
802   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_vec_interp_variant","Variant of algorithm 1-3","HYPRE_BoomerAMGSetInterpVecVariant",jac->vec_interp_variant, &jac->vec_interp_variant,&flg));
803   if (flg) {
804     PetscCallExternal(HYPRE_BoomerAMGSetInterpVecVariant,jac->hsolver,jac->vec_interp_variant);
805   }
806   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_vec_interp_qmax","Max elements per row for each Q","HYPRE_BoomerAMGSetInterpVecQMax",jac->vec_interp_qmax, &jac->vec_interp_qmax,&flg));
807   if (flg) {
808     PetscCallExternal(HYPRE_BoomerAMGSetInterpVecQMax,jac->hsolver,jac->vec_interp_qmax);
809   }
810   PetscCall(PetscOptionsBool("-pc_hypre_boomeramg_vec_interp_smooth","Whether to smooth the interpolation vectors","HYPRE_BoomerAMGSetSmoothInterpVectors",jac->vec_interp_smooth, &jac->vec_interp_smooth,&flg));
811   if (flg) {
812     PetscCallExternal(HYPRE_BoomerAMGSetSmoothInterpVectors,jac->hsolver,jac->vec_interp_smooth);
813   }
814   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_interp_refine","Preprocess the interpolation matrix through iterative weight refinement","HYPRE_BoomerAMGSetInterpRefine",jac->interp_refine, &jac->interp_refine,&flg));
815   if (flg) {
816     PetscCallExternal(HYPRE_BoomerAMGSetInterpRefine,jac->hsolver,jac->interp_refine);
817   }
818   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx,&flg));
819   if (flg) {
820     PetscCallExternal(HYPRE_BoomerAMGSetCycleNumSweeps,jac->hsolver,indx, 1);
821     jac->gridsweeps[0] = indx;
822   }
823   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg));
824   if (flg) {
825     PetscCallExternal(HYPRE_BoomerAMGSetCycleNumSweeps,jac->hsolver,indx, 2);
826     jac->gridsweeps[1] = indx;
827   }
828   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg));
829   if (flg) {
830     PetscCallExternal(HYPRE_BoomerAMGSetCycleNumSweeps,jac->hsolver,indx, 3);
831     jac->gridsweeps[2] = indx;
832   }
833 
834   /* Smooth type */
835   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_smooth_type","Enable more complex smoothers","None",HYPREBoomerAMGSmoothType,PETSC_STATIC_ARRAY_LENGTH(HYPREBoomerAMGSmoothType),HYPREBoomerAMGSmoothType[0],&indx,&flg));
836   if (flg) {
837     jac->smoothtype = indx;
838     PetscCallExternal(HYPRE_BoomerAMGSetSmoothType,jac->hsolver,indx+6);
839     jac->smoothnumlevels = 25;
840     PetscCallExternal(HYPRE_BoomerAMGSetSmoothNumLevels,jac->hsolver,25);
841   }
842 
843   /* Number of smoothing levels */
844   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_smooth_num_levels","Number of levels on which more complex smoothers are used","None",25,&indx,&flg));
845   if (flg && (jac->smoothtype != -1)) {
846     jac->smoothnumlevels = indx;
847     PetscCallExternal(HYPRE_BoomerAMGSetSmoothNumLevels,jac->hsolver,indx);
848   }
849 
850   /* Number of levels for ILU(k) for Euclid */
851   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_eu_level","Number of levels for ILU(k) in Euclid smoother","None",0,&indx,&flg));
852   if (flg && (jac->smoothtype == 3)) {
853     jac->eu_level = indx;
854     PetscCallExternal(HYPRE_BoomerAMGSetEuLevel,jac->hsolver,indx);
855   }
856 
857   /* Filter for ILU(k) for Euclid */
858   double droptolerance;
859   PetscCall(PetscOptionsReal("-pc_hypre_boomeramg_eu_droptolerance","Drop tolerance for ILU(k) in Euclid smoother","None",0,&droptolerance,&flg));
860   if (flg && (jac->smoothtype == 3)) {
861     jac->eu_droptolerance = droptolerance;
862     PetscCallExternal(HYPRE_BoomerAMGSetEuLevel,jac->hsolver,droptolerance);
863   }
864 
865   /* Use Block Jacobi ILUT for Euclid */
866   PetscCall(PetscOptionsBool("-pc_hypre_boomeramg_eu_bj", "Use Block Jacobi for ILU in Euclid smoother?", "None", PETSC_FALSE, &tmp_truth, &flg));
867   if (flg && (jac->smoothtype == 3)) {
868     jac->eu_bj = tmp_truth;
869     PetscCallExternal(HYPRE_BoomerAMGSetEuBJ,jac->hsolver,jac->eu_bj);
870   }
871 
872   /* Relax type */
873   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_relax_type_all","Relax type for the up and down cycles","None",HYPREBoomerAMGRelaxType,PETSC_STATIC_ARRAY_LENGTH(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg));
874   if (flg) {
875     jac->relaxtype[0] = jac->relaxtype[1]  = indx;
876     PetscCallExternal(HYPRE_BoomerAMGSetRelaxType,jac->hsolver, indx);
877     /* by default, coarse type set to 9 */
878     jac->relaxtype[2] = 9;
879     PetscCallExternal(HYPRE_BoomerAMGSetCycleRelaxType,jac->hsolver, 9, 3);
880   }
881   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,PETSC_STATIC_ARRAY_LENGTH(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg));
882   if (flg) {
883     jac->relaxtype[0] = indx;
884     PetscCallExternal(HYPRE_BoomerAMGSetCycleRelaxType,jac->hsolver, indx, 1);
885   }
886   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,PETSC_STATIC_ARRAY_LENGTH(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg));
887   if (flg) {
888     jac->relaxtype[1] = indx;
889     PetscCallExternal(HYPRE_BoomerAMGSetCycleRelaxType,jac->hsolver, indx, 2);
890   }
891   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,PETSC_STATIC_ARRAY_LENGTH(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[9],&indx,&flg));
892   if (flg) {
893     jac->relaxtype[2] = indx;
894     PetscCallExternal(HYPRE_BoomerAMGSetCycleRelaxType,jac->hsolver, indx, 3);
895   }
896 
897   /* Relaxation Weight */
898   PetscCall(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));
899   if (flg) {
900     PetscCallExternal(HYPRE_BoomerAMGSetRelaxWt,jac->hsolver,tmpdbl);
901     jac->relaxweight = tmpdbl;
902   }
903 
904   n         = 2;
905   twodbl[0] = twodbl[1] = 1.0;
906   PetscCall(PetscOptionsRealArray("-pc_hypre_boomeramg_relax_weight_level","Set the relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg));
907   if (flg) {
908     if (n == 2) {
909       indx =  (int)PetscAbsReal(twodbl[1]);
910       PetscCallExternal(HYPRE_BoomerAMGSetLevelRelaxWt,jac->hsolver,twodbl[0],indx);
911     } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relax weight level: you must provide 2 values separated by a comma (and no space), you provided %" PetscInt_FMT,n);
912   }
913 
914   /* Outer relaxation Weight */
915   PetscCall(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));
916   if (flg) {
917     PetscCallExternal(HYPRE_BoomerAMGSetOuterWt,jac->hsolver, tmpdbl);
918     jac->outerrelaxweight = tmpdbl;
919   }
920 
921   n         = 2;
922   twodbl[0] = twodbl[1] = 1.0;
923   PetscCall(PetscOptionsRealArray("-pc_hypre_boomeramg_outer_relax_weight_level","Set the outer relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg));
924   if (flg) {
925     if (n == 2) {
926       indx =  (int)PetscAbsReal(twodbl[1]);
927       PetscCallExternal(HYPRE_BoomerAMGSetLevelOuterWt,jac->hsolver, twodbl[0], indx);
928     } else SETERRQ(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 %" PetscInt_FMT,n);
929   }
930 
931   /* the Relax Order */
932   PetscCall(PetscOptionsBool("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg));
933 
934   if (flg && tmp_truth) {
935     jac->relaxorder = 0;
936     PetscCallExternal(HYPRE_BoomerAMGSetRelaxOrder,jac->hsolver, jac->relaxorder);
937   }
938   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,PETSC_STATIC_ARRAY_LENGTH(HYPREBoomerAMGMeasureType),HYPREBoomerAMGMeasureType[0],&indx,&flg));
939   if (flg) {
940     jac->measuretype = indx;
941     PetscCallExternal(HYPRE_BoomerAMGSetMeasureType,jac->hsolver,jac->measuretype);
942   }
943   /* update list length 3/07 */
944   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,PETSC_STATIC_ARRAY_LENGTH(HYPREBoomerAMGCoarsenType),HYPREBoomerAMGCoarsenType[6],&indx,&flg));
945   if (flg) {
946     jac->coarsentype = indx;
947     PetscCallExternal(HYPRE_BoomerAMGSetCoarsenType,jac->hsolver,jac->coarsentype);
948   }
949 
950   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_max_coarse_size", "Maximum size of coarsest grid", "None", jac->maxc, &jac->maxc, &flg));
951   if (flg) {
952     PetscCallExternal(HYPRE_BoomerAMGSetMaxCoarseSize,jac->hsolver, jac->maxc);
953   }
954   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_min_coarse_size", "Minimum size of coarsest grid", "None", jac->minc, &jac->minc, &flg));
955   if (flg) {
956     PetscCallExternal(HYPRE_BoomerAMGSetMinCoarseSize,jac->hsolver, jac->minc);
957   }
958 #if PETSC_PKG_HYPRE_VERSION_GE(2,23,0)
959   // global parameter but is closely associated with BoomerAMG
960   PetscCall(PetscOptionsEList("-pc_mg_galerkin_mat_product_algorithm","Type of SpGEMM to use in hypre (only for now)","PCMGGalerkinSetMatProductAlgorithm",PCHYPRESpgemmTypes,PETSC_STATIC_ARRAY_LENGTH(PCHYPRESpgemmTypes),PCHYPRESpgemmTypes[0],&indx,&flg));
961   if (!flg) indx = 0;
962   PetscCall(PCMGGalerkinSetMatProductAlgorithm_HYPRE_BoomerAMG(pc,PCHYPRESpgemmTypes[indx]));
963 #endif
964   /* AIR */
965 #if PETSC_PKG_HYPRE_VERSION_GE(2,18,0)
966   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_restriction_type", "Type of AIR method (distance 1 or 2, 0 means no AIR)", "None", jac->Rtype, &jac->Rtype, NULL));
967   PetscCallExternal(HYPRE_BoomerAMGSetRestriction,jac->hsolver,jac->Rtype);
968   if (jac->Rtype) {
969     jac->interptype = 100; /* no way we can pass this with strings... Set it as default as in MFEM, then users can still customize it back to a different one */
970 
971     PetscCall(PetscOptionsReal("-pc_hypre_boomeramg_strongthresholdR","Threshold for R","None",jac->Rstrongthreshold,&jac->Rstrongthreshold,NULL));
972     PetscCallExternal(HYPRE_BoomerAMGSetStrongThresholdR,jac->hsolver,jac->Rstrongthreshold);
973 
974     PetscCall(PetscOptionsReal("-pc_hypre_boomeramg_filterthresholdR","Filter threshold for R","None",jac->Rfilterthreshold,&jac->Rfilterthreshold,NULL));
975     PetscCallExternal(HYPRE_BoomerAMGSetFilterThresholdR,jac->hsolver,jac->Rfilterthreshold);
976 
977     PetscCall(PetscOptionsReal("-pc_hypre_boomeramg_Adroptol","Defines the drop tolerance for the A-matrices from the 2nd level of AMG","None",jac->Adroptol,&jac->Adroptol,NULL));
978     PetscCallExternal(HYPRE_BoomerAMGSetADropTol,jac->hsolver,jac->Adroptol);
979 
980     PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_Adroptype","Drops the entries that are not on the diagonal and smaller than its row norm: type 1: 1-norm, 2: 2-norm, -1: infinity norm","None",jac->Adroptype,&jac->Adroptype,NULL));
981     PetscCallExternal(HYPRE_BoomerAMGSetADropType,jac->hsolver,jac->Adroptype);
982   }
983 #endif
984 
985 #if PETSC_PKG_HYPRE_VERSION_LE(9,9,9)
986   PetscCheck(!jac->Rtype || !jac->agg_nl,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"-pc_hypre_boomeramg_restriction_type (%" PetscInt_FMT ") and -pc_hypre_boomeramg_agg_nl (%" PetscInt_FMT ")",jac->Rtype,jac->agg_nl);
987 #endif
988 
989   /* new 3/07 */
990   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,PETSC_STATIC_ARRAY_LENGTH(HYPREBoomerAMGInterpType),HYPREBoomerAMGInterpType[0],&indx,&flg));
991   if (flg || jac->Rtype) {
992     if (flg) jac->interptype = indx;
993     PetscCallExternal(HYPRE_BoomerAMGSetInterpType,jac->hsolver,jac->interptype);
994   }
995 
996   PetscCall(PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg));
997   if (flg) {
998     level = 3;
999     PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,NULL));
1000 
1001     jac->printstatistics = PETSC_TRUE;
1002     PetscCallExternal(HYPRE_BoomerAMGSetPrintLevel,jac->hsolver,level);
1003   }
1004 
1005   PetscCall(PetscOptionsName("-pc_hypre_boomeramg_print_debug","Print debug information","None",&flg));
1006   if (flg) {
1007     level = 3;
1008     PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,NULL));
1009 
1010     jac->printstatistics = PETSC_TRUE;
1011     PetscCallExternal(HYPRE_BoomerAMGSetDebugFlag,jac->hsolver,level);
1012   }
1013 
1014   PetscCall(PetscOptionsBool("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg));
1015   if (flg && tmp_truth) {
1016     PetscInt tmp_int;
1017     PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg));
1018     if (flg) jac->nodal_relax_levels = tmp_int;
1019     PetscCallExternal(HYPRE_BoomerAMGSetSmoothType,jac->hsolver,6);
1020     PetscCallExternal(HYPRE_BoomerAMGSetDomainType,jac->hsolver,1);
1021     PetscCallExternal(HYPRE_BoomerAMGSetOverlap,jac->hsolver,0);
1022     PetscCallExternal(HYPRE_BoomerAMGSetSmoothNumLevels,jac->hsolver,jac->nodal_relax_levels);
1023   }
1024 
1025   PetscCall(PetscOptionsBool("-pc_hypre_boomeramg_keeptranspose", "Avoid transpose matvecs in preconditioner application", "None", jac->keeptranspose, &jac->keeptranspose, NULL));
1026   PetscCallExternal(HYPRE_BoomerAMGSetKeepTranspose,jac->hsolver,jac->keeptranspose ? 1 : 0);
1027 
1028   /* options for ParaSails solvers */
1029   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,PETSC_STATIC_ARRAY_LENGTH(symtlist),symtlist[0],&indx,&flg));
1030   if (flg) {
1031     jac->symt = indx;
1032     PetscCallExternal(HYPRE_BoomerAMGSetSym,jac->hsolver,jac->symt);
1033   }
1034 
1035   PetscOptionsHeadEnd();
1036   PetscFunctionReturn(0);
1037 }
1038 
1039 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)
1040 {
1041   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1042   HYPRE_Int      oits;
1043 
1044   PetscFunctionBegin;
1045   PetscCall(PetscCitationsRegister(hypreCitation,&cite));
1046   PetscCallExternal(HYPRE_BoomerAMGSetMaxIter,jac->hsolver,its*jac->maxiter);
1047   PetscCallExternal(HYPRE_BoomerAMGSetTol,jac->hsolver,rtol);
1048   jac->applyrichardson = PETSC_TRUE;
1049   PetscCall(PCApply_HYPRE(pc,b,y));
1050   jac->applyrichardson = PETSC_FALSE;
1051   PetscCallExternal(HYPRE_BoomerAMGGetNumIterations,jac->hsolver,&oits);
1052   *outits = oits;
1053   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
1054   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
1055   PetscCallExternal(HYPRE_BoomerAMGSetTol,jac->hsolver,jac->tol);
1056   PetscCallExternal(HYPRE_BoomerAMGSetMaxIter,jac->hsolver,jac->maxiter);
1057   PetscFunctionReturn(0);
1058 }
1059 
1060 static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
1061 {
1062   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1063   PetscBool      iascii;
1064 
1065   PetscFunctionBegin;
1066   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
1067   if (iascii) {
1068     PetscCall(PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n"));
1069     PetscCall(PetscViewerASCIIPrintf(viewer,"    Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]));
1070     PetscCall(PetscViewerASCIIPrintf(viewer,"    Maximum number of levels %" PetscInt_FMT "\n",jac->maxlevels));
1071     PetscCall(PetscViewerASCIIPrintf(viewer,"    Maximum number of iterations PER hypre call %" PetscInt_FMT "\n",jac->maxiter));
1072     PetscCall(PetscViewerASCIIPrintf(viewer,"    Convergence tolerance PER hypre call %g\n",(double)jac->tol));
1073     PetscCall(PetscViewerASCIIPrintf(viewer,"    Threshold for strong coupling %g\n",(double)jac->strongthreshold));
1074     PetscCall(PetscViewerASCIIPrintf(viewer,"    Interpolation truncation factor %g\n",(double)jac->truncfactor));
1075     PetscCall(PetscViewerASCIIPrintf(viewer,"    Interpolation: max elements per row %" PetscInt_FMT "\n",jac->pmax));
1076     if (jac->interp_refine) {
1077       PetscCall(PetscViewerASCIIPrintf(viewer,"    Interpolation: number of steps of weighted refinement %" PetscInt_FMT "\n",jac->interp_refine));
1078     }
1079     PetscCall(PetscViewerASCIIPrintf(viewer,"    Number of levels of aggressive coarsening %" PetscInt_FMT "\n",jac->agg_nl));
1080     PetscCall(PetscViewerASCIIPrintf(viewer,"    Number of paths for aggressive coarsening %" PetscInt_FMT "\n",jac->agg_num_paths));
1081 
1082     PetscCall(PetscViewerASCIIPrintf(viewer,"    Maximum row sums %g\n",(double)jac->maxrowsum));
1083 
1084     PetscCall(PetscViewerASCIIPrintf(viewer,"    Sweeps down         %" PetscInt_FMT "\n",jac->gridsweeps[0]));
1085     PetscCall(PetscViewerASCIIPrintf(viewer,"    Sweeps up           %" PetscInt_FMT "\n",jac->gridsweeps[1]));
1086     PetscCall(PetscViewerASCIIPrintf(viewer,"    Sweeps on coarse    %" PetscInt_FMT "\n",jac->gridsweeps[2]));
1087 
1088     PetscCall(PetscViewerASCIIPrintf(viewer,"    Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]));
1089     PetscCall(PetscViewerASCIIPrintf(viewer,"    Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]));
1090     PetscCall(PetscViewerASCIIPrintf(viewer,"    Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]));
1091 
1092     PetscCall(PetscViewerASCIIPrintf(viewer,"    Relax weight  (all)      %g\n",(double)jac->relaxweight));
1093     PetscCall(PetscViewerASCIIPrintf(viewer,"    Outer relax weight (all) %g\n",(double)jac->outerrelaxweight));
1094 
1095     if (jac->relaxorder) {
1096       PetscCall(PetscViewerASCIIPrintf(viewer,"    Using CF-relaxation\n"));
1097     } else {
1098       PetscCall(PetscViewerASCIIPrintf(viewer,"    Not using CF-relaxation\n"));
1099     }
1100     if (jac->smoothtype!=-1) {
1101       PetscCall(PetscViewerASCIIPrintf(viewer,"    Smooth type          %s\n",HYPREBoomerAMGSmoothType[jac->smoothtype]));
1102       PetscCall(PetscViewerASCIIPrintf(viewer,"    Smooth num levels    %" PetscInt_FMT "\n",jac->smoothnumlevels));
1103     } else {
1104       PetscCall(PetscViewerASCIIPrintf(viewer,"    Not using more complex smoothers.\n"));
1105     }
1106     if (jac->smoothtype==3) {
1107       PetscCall(PetscViewerASCIIPrintf(viewer,"    Euclid ILU(k) levels %" PetscInt_FMT "\n",jac->eu_level));
1108       PetscCall(PetscViewerASCIIPrintf(viewer,"    Euclid ILU(k) drop tolerance %g\n",(double)jac->eu_droptolerance));
1109       PetscCall(PetscViewerASCIIPrintf(viewer,"    Euclid ILU use Block-Jacobi? %" PetscInt_FMT "\n",jac->eu_bj));
1110     }
1111     PetscCall(PetscViewerASCIIPrintf(viewer,"    Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]));
1112     PetscCall(PetscViewerASCIIPrintf(viewer,"    Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]));
1113     PetscCall(PetscViewerASCIIPrintf(viewer,"    Interpolation type  %s\n",jac->interptype != 100 ? HYPREBoomerAMGInterpType[jac->interptype] : "1pt"));
1114     if (jac->nodal_coarsening) {
1115       PetscCall(PetscViewerASCIIPrintf(viewer,"    Using nodal coarsening with HYPRE_BOOMERAMGSetNodal() %" PetscInt_FMT "\n",jac->nodal_coarsening));
1116     }
1117     if (jac->vec_interp_variant) {
1118       PetscCall(PetscViewerASCIIPrintf(viewer,"    HYPRE_BoomerAMGSetInterpVecVariant() %" PetscInt_FMT "\n",jac->vec_interp_variant));
1119       PetscCall(PetscViewerASCIIPrintf(viewer,"    HYPRE_BoomerAMGSetInterpVecQMax() %" PetscInt_FMT "\n",jac->vec_interp_qmax));
1120       PetscCall(PetscViewerASCIIPrintf(viewer,"    HYPRE_BoomerAMGSetSmoothInterpVectors() %d\n",jac->vec_interp_smooth));
1121     }
1122     if (jac->nodal_relax) {
1123       PetscCall(PetscViewerASCIIPrintf(viewer,"    Using nodal relaxation via Schwarz smoothing on levels %" PetscInt_FMT "\n",jac->nodal_relax_levels));
1124     }
1125 #if PETSC_PKG_HYPRE_VERSION_GE(2,23,0)
1126     PetscCall(PetscViewerASCIIPrintf(viewer,"    SpGEMM type         %s\n",jac->spgemm_type));
1127 #endif
1128     /* AIR */
1129     if (jac->Rtype) {
1130       PetscCall(PetscViewerASCIIPrintf(viewer,"    Using approximate ideal restriction type %" PetscInt_FMT "\n",jac->Rtype));
1131       PetscCall(PetscViewerASCIIPrintf(viewer,"      Threshold for R %g\n",(double)jac->Rstrongthreshold));
1132       PetscCall(PetscViewerASCIIPrintf(viewer,"      Filter for R %g\n",(double)jac->Rfilterthreshold));
1133       PetscCall(PetscViewerASCIIPrintf(viewer,"      A drop tolerance %g\n",(double)jac->Adroptol));
1134       PetscCall(PetscViewerASCIIPrintf(viewer,"      A drop type %" PetscInt_FMT "\n",jac->Adroptype));
1135     }
1136   }
1137   PetscFunctionReturn(0);
1138 }
1139 
1140 /* --------------------------------------------------------------------------------------------*/
1141 static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PetscOptionItems *PetscOptionsObject,PC pc)
1142 {
1143   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1144   PetscInt       indx;
1145   PetscBool      flag;
1146   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};
1147 
1148   PetscFunctionBegin;
1149   PetscOptionsHeadBegin(PetscOptionsObject,"HYPRE ParaSails Options");
1150   PetscCall(PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0));
1151   PetscCall(PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshold,&jac->threshold,&flag));
1152   if (flag) PetscCallExternal(HYPRE_ParaSailsSetParams,jac->hsolver,jac->threshold,jac->nlevels);
1153 
1154   PetscCall(PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag));
1155   if (flag) PetscCallExternal(HYPRE_ParaSailsSetFilter,jac->hsolver,jac->filter);
1156 
1157   PetscCall(PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag));
1158   if (flag) PetscCallExternal(HYPRE_ParaSailsSetLoadbal,jac->hsolver,jac->loadbal);
1159 
1160   PetscCall(PetscOptionsBool("-pc_hypre_parasails_logging","Print info to screen","None",(PetscBool)jac->logging,(PetscBool*)&jac->logging,&flag));
1161   if (flag) PetscCallExternal(HYPRE_ParaSailsSetLogging,jac->hsolver,jac->logging);
1162 
1163   PetscCall(PetscOptionsBool("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscBool)jac->ruse,(PetscBool*)&jac->ruse,&flag));
1164   if (flag) PetscCallExternal(HYPRE_ParaSailsSetReuse,jac->hsolver,jac->ruse);
1165 
1166   PetscCall(PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,PETSC_STATIC_ARRAY_LENGTH(symtlist),symtlist[0],&indx,&flag));
1167   if (flag) {
1168     jac->symt = indx;
1169     PetscCallExternal(HYPRE_ParaSailsSetSym,jac->hsolver,jac->symt);
1170   }
1171 
1172   PetscOptionsHeadEnd();
1173   PetscFunctionReturn(0);
1174 }
1175 
1176 static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
1177 {
1178   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1179   PetscBool      iascii;
1180   const char     *symt = 0;
1181 
1182   PetscFunctionBegin;
1183   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
1184   if (iascii) {
1185     PetscCall(PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n"));
1186     PetscCall(PetscViewerASCIIPrintf(viewer,"    nlevels %" PetscInt_FMT "\n",jac->nlevels));
1187     PetscCall(PetscViewerASCIIPrintf(viewer,"    threshold %g\n",(double)jac->threshold));
1188     PetscCall(PetscViewerASCIIPrintf(viewer,"    filter %g\n",(double)jac->filter));
1189     PetscCall(PetscViewerASCIIPrintf(viewer,"    load balance %g\n",(double)jac->loadbal));
1190     PetscCall(PetscViewerASCIIPrintf(viewer,"    reuse nonzero structure %s\n",PetscBools[jac->ruse]));
1191     PetscCall(PetscViewerASCIIPrintf(viewer,"    print info to screen %s\n",PetscBools[jac->logging]));
1192     if (!jac->symt) symt = "nonsymmetric matrix and preconditioner";
1193     else if (jac->symt == 1) symt = "SPD matrix and preconditioner";
1194     else if (jac->symt == 2) symt = "nonsymmetric matrix but SPD preconditioner";
1195     else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %" PetscInt_FMT,jac->symt);
1196     PetscCall(PetscViewerASCIIPrintf(viewer,"    %s\n",symt));
1197   }
1198   PetscFunctionReturn(0);
1199 }
1200 /* --------------------------------------------------------------------------------------------*/
1201 static PetscErrorCode PCSetFromOptions_HYPRE_AMS(PetscOptionItems *PetscOptionsObject,PC pc)
1202 {
1203   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1204   PetscInt       n;
1205   PetscBool      flag,flag2,flag3,flag4;
1206 
1207   PetscFunctionBegin;
1208   PetscOptionsHeadBegin(PetscOptionsObject,"HYPRE AMS Options");
1209   PetscCall(PetscOptionsInt("-pc_hypre_ams_print_level","Debugging output level for AMS","None",jac->as_print,&jac->as_print,&flag));
1210   if (flag) PetscCallExternal(HYPRE_AMSSetPrintLevel,jac->hsolver,jac->as_print);
1211   PetscCall(PetscOptionsInt("-pc_hypre_ams_max_iter","Maximum number of AMS multigrid iterations within PCApply","None",jac->as_max_iter,&jac->as_max_iter,&flag));
1212   if (flag) PetscCallExternal(HYPRE_AMSSetMaxIter,jac->hsolver,jac->as_max_iter);
1213   PetscCall(PetscOptionsInt("-pc_hypre_ams_cycle_type","Cycle type for AMS multigrid","None",jac->ams_cycle_type,&jac->ams_cycle_type,&flag));
1214   if (flag) PetscCallExternal(HYPRE_AMSSetCycleType,jac->hsolver,jac->ams_cycle_type);
1215   PetscCall(PetscOptionsReal("-pc_hypre_ams_tol","Error tolerance for AMS multigrid","None",jac->as_tol,&jac->as_tol,&flag));
1216   if (flag) PetscCallExternal(HYPRE_AMSSetTol,jac->hsolver,jac->as_tol);
1217   PetscCall(PetscOptionsInt("-pc_hypre_ams_relax_type","Relaxation type for AMS smoother","None",jac->as_relax_type,&jac->as_relax_type,&flag));
1218   PetscCall(PetscOptionsInt("-pc_hypre_ams_relax_times","Number of relaxation steps for AMS smoother","None",jac->as_relax_times,&jac->as_relax_times,&flag2));
1219   PetscCall(PetscOptionsReal("-pc_hypre_ams_relax_weight","Relaxation weight for AMS smoother","None",jac->as_relax_weight,&jac->as_relax_weight,&flag3));
1220   PetscCall(PetscOptionsReal("-pc_hypre_ams_omega","SSOR coefficient for AMS smoother","None",jac->as_omega,&jac->as_omega,&flag4));
1221   if (flag || flag2 || flag3 || flag4) {
1222     PetscCallExternal(HYPRE_AMSSetSmoothingOptions,jac->hsolver,jac->as_relax_type,
1223                                                                       jac->as_relax_times,
1224                                                                       jac->as_relax_weight,
1225                                                                       jac->as_omega);
1226   }
1227   PetscCall(PetscOptionsReal("-pc_hypre_ams_amg_alpha_theta","Threshold for strong coupling of vector Poisson AMG solver","None",jac->as_amg_alpha_theta,&jac->as_amg_alpha_theta,&flag));
1228   n = 5;
1229   PetscCall(PetscOptionsIntArray("-pc_hypre_ams_amg_alpha_options","AMG options for vector Poisson","None",jac->as_amg_alpha_opts,&n,&flag2));
1230   if (flag || flag2) {
1231     PetscCallExternal(HYPRE_AMSSetAlphaAMGOptions,jac->hsolver,jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
1232                                                                      jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
1233                                                                      jac->as_amg_alpha_opts[2],       /* AMG relax_type */
1234                                                                      jac->as_amg_alpha_theta,
1235                                                                      jac->as_amg_alpha_opts[3],       /* AMG interp_type */
1236                                                                      jac->as_amg_alpha_opts[4]);     /* AMG Pmax */
1237   }
1238   PetscCall(PetscOptionsReal("-pc_hypre_ams_amg_beta_theta","Threshold for strong coupling of scalar Poisson AMG solver","None",jac->as_amg_beta_theta,&jac->as_amg_beta_theta,&flag));
1239   n = 5;
1240   PetscCall(PetscOptionsIntArray("-pc_hypre_ams_amg_beta_options","AMG options for scalar Poisson solver","None",jac->as_amg_beta_opts,&n,&flag2));
1241   if (flag || flag2) {
1242     PetscCallExternal(HYPRE_AMSSetBetaAMGOptions,jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
1243                                                                     jac->as_amg_beta_opts[1],       /* AMG agg_levels */
1244                                                                     jac->as_amg_beta_opts[2],       /* AMG relax_type */
1245                                                                     jac->as_amg_beta_theta,
1246                                                                     jac->as_amg_beta_opts[3],       /* AMG interp_type */
1247                                                                     jac->as_amg_beta_opts[4]);     /* AMG Pmax */
1248   }
1249   PetscCall(PetscOptionsInt("-pc_hypre_ams_projection_frequency","Frequency at which a projection onto the compatible subspace for problems with zero conductivity regions is performed","None",jac->ams_proj_freq,&jac->ams_proj_freq,&flag));
1250   if (flag) { /* override HYPRE's default only if the options is used */
1251     PetscCallExternal(HYPRE_AMSSetProjectionFrequency,jac->hsolver,jac->ams_proj_freq);
1252   }
1253   PetscOptionsHeadEnd();
1254   PetscFunctionReturn(0);
1255 }
1256 
1257 static PetscErrorCode PCView_HYPRE_AMS(PC pc,PetscViewer viewer)
1258 {
1259   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1260   PetscBool      iascii;
1261 
1262   PetscFunctionBegin;
1263   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
1264   if (iascii) {
1265     PetscCall(PetscViewerASCIIPrintf(viewer,"  HYPRE AMS preconditioning\n"));
1266     PetscCall(PetscViewerASCIIPrintf(viewer,"    subspace iterations per application %" PetscInt_FMT "\n",jac->as_max_iter));
1267     PetscCall(PetscViewerASCIIPrintf(viewer,"    subspace cycle type %" PetscInt_FMT "\n",jac->ams_cycle_type));
1268     PetscCall(PetscViewerASCIIPrintf(viewer,"    subspace iteration tolerance %g\n",(double)jac->as_tol));
1269     PetscCall(PetscViewerASCIIPrintf(viewer,"    smoother type %" PetscInt_FMT "\n",jac->as_relax_type));
1270     PetscCall(PetscViewerASCIIPrintf(viewer,"    number of smoothing steps %" PetscInt_FMT "\n",jac->as_relax_times));
1271     PetscCall(PetscViewerASCIIPrintf(viewer,"    smoother weight %g\n",(double)jac->as_relax_weight));
1272     PetscCall(PetscViewerASCIIPrintf(viewer,"    smoother omega %g\n",(double)jac->as_omega));
1273     if (jac->alpha_Poisson) {
1274       PetscCall(PetscViewerASCIIPrintf(viewer,"    vector Poisson solver (passed in by user)\n"));
1275     } else {
1276       PetscCall(PetscViewerASCIIPrintf(viewer,"    vector Poisson solver (computed) \n"));
1277     }
1278     PetscCall(PetscViewerASCIIPrintf(viewer,"        boomerAMG coarsening type %" PetscInt_FMT "\n",jac->as_amg_alpha_opts[0]));
1279     PetscCall(PetscViewerASCIIPrintf(viewer,"        boomerAMG levels of aggressive coarsening %" PetscInt_FMT "\n",jac->as_amg_alpha_opts[1]));
1280     PetscCall(PetscViewerASCIIPrintf(viewer,"        boomerAMG relaxation type %" PetscInt_FMT "\n",jac->as_amg_alpha_opts[2]));
1281     PetscCall(PetscViewerASCIIPrintf(viewer,"        boomerAMG interpolation type %" PetscInt_FMT "\n",jac->as_amg_alpha_opts[3]));
1282     PetscCall(PetscViewerASCIIPrintf(viewer,"        boomerAMG max nonzero elements in interpolation rows %" PetscInt_FMT "\n",jac->as_amg_alpha_opts[4]));
1283     PetscCall(PetscViewerASCIIPrintf(viewer,"        boomerAMG strength threshold %g\n",(double)jac->as_amg_alpha_theta));
1284     if (!jac->ams_beta_is_zero) {
1285       if (jac->beta_Poisson) {
1286         PetscCall(PetscViewerASCIIPrintf(viewer,"    scalar Poisson solver (passed in by user)\n"));
1287       } else {
1288         PetscCall(PetscViewerASCIIPrintf(viewer,"    scalar Poisson solver (computed) \n"));
1289       }
1290       PetscCall(PetscViewerASCIIPrintf(viewer,"        boomerAMG coarsening type %" PetscInt_FMT "\n",jac->as_amg_beta_opts[0]));
1291       PetscCall(PetscViewerASCIIPrintf(viewer,"        boomerAMG levels of aggressive coarsening %" PetscInt_FMT "\n",jac->as_amg_beta_opts[1]));
1292       PetscCall(PetscViewerASCIIPrintf(viewer,"        boomerAMG relaxation type %" PetscInt_FMT "\n",jac->as_amg_beta_opts[2]));
1293       PetscCall(PetscViewerASCIIPrintf(viewer,"        boomerAMG interpolation type %" PetscInt_FMT "\n",jac->as_amg_beta_opts[3]));
1294       PetscCall(PetscViewerASCIIPrintf(viewer,"        boomerAMG max nonzero elements in interpolation rows %" PetscInt_FMT "\n",jac->as_amg_beta_opts[4]));
1295       PetscCall(PetscViewerASCIIPrintf(viewer,"        boomerAMG strength threshold %g\n",(double)jac->as_amg_beta_theta));
1296       if (jac->ams_beta_is_zero_part) {
1297         PetscCall(PetscViewerASCIIPrintf(viewer,"        compatible subspace projection frequency %" PetscInt_FMT " (-1 HYPRE uses default)\n",jac->ams_proj_freq));
1298       }
1299     } else {
1300       PetscCall(PetscViewerASCIIPrintf(viewer,"    scalar Poisson solver not used (zero-conductivity everywhere) \n"));
1301     }
1302   }
1303   PetscFunctionReturn(0);
1304 }
1305 
1306 static PetscErrorCode PCSetFromOptions_HYPRE_ADS(PetscOptionItems *PetscOptionsObject,PC pc)
1307 {
1308   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1309   PetscInt       n;
1310   PetscBool      flag,flag2,flag3,flag4;
1311 
1312   PetscFunctionBegin;
1313   PetscOptionsHeadBegin(PetscOptionsObject,"HYPRE ADS Options");
1314   PetscCall(PetscOptionsInt("-pc_hypre_ads_print_level","Debugging output level for ADS","None",jac->as_print,&jac->as_print,&flag));
1315   if (flag) PetscCallExternal(HYPRE_ADSSetPrintLevel,jac->hsolver,jac->as_print);
1316   PetscCall(PetscOptionsInt("-pc_hypre_ads_max_iter","Maximum number of ADS multigrid iterations within PCApply","None",jac->as_max_iter,&jac->as_max_iter,&flag));
1317   if (flag) PetscCallExternal(HYPRE_ADSSetMaxIter,jac->hsolver,jac->as_max_iter);
1318   PetscCall(PetscOptionsInt("-pc_hypre_ads_cycle_type","Cycle type for ADS multigrid","None",jac->ads_cycle_type,&jac->ads_cycle_type,&flag));
1319   if (flag) PetscCallExternal(HYPRE_ADSSetCycleType,jac->hsolver,jac->ads_cycle_type);
1320   PetscCall(PetscOptionsReal("-pc_hypre_ads_tol","Error tolerance for ADS multigrid","None",jac->as_tol,&jac->as_tol,&flag));
1321   if (flag) PetscCallExternal(HYPRE_ADSSetTol,jac->hsolver,jac->as_tol);
1322   PetscCall(PetscOptionsInt("-pc_hypre_ads_relax_type","Relaxation type for ADS smoother","None",jac->as_relax_type,&jac->as_relax_type,&flag));
1323   PetscCall(PetscOptionsInt("-pc_hypre_ads_relax_times","Number of relaxation steps for ADS smoother","None",jac->as_relax_times,&jac->as_relax_times,&flag2));
1324   PetscCall(PetscOptionsReal("-pc_hypre_ads_relax_weight","Relaxation weight for ADS smoother","None",jac->as_relax_weight,&jac->as_relax_weight,&flag3));
1325   PetscCall(PetscOptionsReal("-pc_hypre_ads_omega","SSOR coefficient for ADS smoother","None",jac->as_omega,&jac->as_omega,&flag4));
1326   if (flag || flag2 || flag3 || flag4) {
1327     PetscCallExternal(HYPRE_ADSSetSmoothingOptions,jac->hsolver,jac->as_relax_type,
1328                                                                       jac->as_relax_times,
1329                                                                       jac->as_relax_weight,
1330                                                                       jac->as_omega);
1331   }
1332   PetscCall(PetscOptionsReal("-pc_hypre_ads_ams_theta","Threshold for strong coupling of AMS solver inside ADS","None",jac->as_amg_alpha_theta,&jac->as_amg_alpha_theta,&flag));
1333   n = 5;
1334   PetscCall(PetscOptionsIntArray("-pc_hypre_ads_ams_options","AMG options for AMS solver inside ADS","None",jac->as_amg_alpha_opts,&n,&flag2));
1335   PetscCall(PetscOptionsInt("-pc_hypre_ads_ams_cycle_type","Cycle type for AMS solver inside ADS","None",jac->ams_cycle_type,&jac->ams_cycle_type,&flag3));
1336   if (flag || flag2 || flag3) {
1337     PetscCallExternal(HYPRE_ADSSetAMSOptions,jac->hsolver,jac->ams_cycle_type,             /* AMS cycle type */
1338                                                                 jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
1339                                                                 jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
1340                                                                 jac->as_amg_alpha_opts[2],       /* AMG relax_type */
1341                                                                 jac->as_amg_alpha_theta,
1342                                                                 jac->as_amg_alpha_opts[3],       /* AMG interp_type */
1343                                                                 jac->as_amg_alpha_opts[4]);     /* AMG Pmax */
1344   }
1345   PetscCall(PetscOptionsReal("-pc_hypre_ads_amg_theta","Threshold for strong coupling of vector AMG solver inside ADS","None",jac->as_amg_beta_theta,&jac->as_amg_beta_theta,&flag));
1346   n = 5;
1347   PetscCall(PetscOptionsIntArray("-pc_hypre_ads_amg_options","AMG options for vector AMG solver inside ADS","None",jac->as_amg_beta_opts,&n,&flag2));
1348   if (flag || flag2) {
1349     PetscCallExternal(HYPRE_ADSSetAMGOptions,jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
1350                                                                 jac->as_amg_beta_opts[1],       /* AMG agg_levels */
1351                                                                 jac->as_amg_beta_opts[2],       /* AMG relax_type */
1352                                                                 jac->as_amg_beta_theta,
1353                                                                 jac->as_amg_beta_opts[3],       /* AMG interp_type */
1354                                                                 jac->as_amg_beta_opts[4]);     /* AMG Pmax */
1355   }
1356   PetscOptionsHeadEnd();
1357   PetscFunctionReturn(0);
1358 }
1359 
1360 static PetscErrorCode PCView_HYPRE_ADS(PC pc,PetscViewer viewer)
1361 {
1362   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1363   PetscBool      iascii;
1364 
1365   PetscFunctionBegin;
1366   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
1367   if (iascii) {
1368     PetscCall(PetscViewerASCIIPrintf(viewer,"  HYPRE ADS preconditioning\n"));
1369     PetscCall(PetscViewerASCIIPrintf(viewer,"    subspace iterations per application %" PetscInt_FMT "\n",jac->as_max_iter));
1370     PetscCall(PetscViewerASCIIPrintf(viewer,"    subspace cycle type %" PetscInt_FMT "\n",jac->ads_cycle_type));
1371     PetscCall(PetscViewerASCIIPrintf(viewer,"    subspace iteration tolerance %g\n",(double)jac->as_tol));
1372     PetscCall(PetscViewerASCIIPrintf(viewer,"    smoother type %" PetscInt_FMT "\n",jac->as_relax_type));
1373     PetscCall(PetscViewerASCIIPrintf(viewer,"    number of smoothing steps %" PetscInt_FMT "\n",jac->as_relax_times));
1374     PetscCall(PetscViewerASCIIPrintf(viewer,"    smoother weight %g\n",(double)jac->as_relax_weight));
1375     PetscCall(PetscViewerASCIIPrintf(viewer,"    smoother omega %g\n",(double)jac->as_omega));
1376     PetscCall(PetscViewerASCIIPrintf(viewer,"    AMS solver using boomerAMG\n"));
1377     PetscCall(PetscViewerASCIIPrintf(viewer,"        subspace cycle type %" PetscInt_FMT "\n",jac->ams_cycle_type));
1378     PetscCall(PetscViewerASCIIPrintf(viewer,"        coarsening type %" PetscInt_FMT "\n",jac->as_amg_alpha_opts[0]));
1379     PetscCall(PetscViewerASCIIPrintf(viewer,"        levels of aggressive coarsening %" PetscInt_FMT "\n",jac->as_amg_alpha_opts[1]));
1380     PetscCall(PetscViewerASCIIPrintf(viewer,"        relaxation type %" PetscInt_FMT "\n",jac->as_amg_alpha_opts[2]));
1381     PetscCall(PetscViewerASCIIPrintf(viewer,"        interpolation type %" PetscInt_FMT "\n",jac->as_amg_alpha_opts[3]));
1382     PetscCall(PetscViewerASCIIPrintf(viewer,"        max nonzero elements in interpolation rows %" PetscInt_FMT "\n",jac->as_amg_alpha_opts[4]));
1383     PetscCall(PetscViewerASCIIPrintf(viewer,"        strength threshold %g\n",(double)jac->as_amg_alpha_theta));
1384     PetscCall(PetscViewerASCIIPrintf(viewer,"    vector Poisson solver using boomerAMG\n"));
1385     PetscCall(PetscViewerASCIIPrintf(viewer,"        coarsening type %" PetscInt_FMT "\n",jac->as_amg_beta_opts[0]));
1386     PetscCall(PetscViewerASCIIPrintf(viewer,"        levels of aggressive coarsening %" PetscInt_FMT "\n",jac->as_amg_beta_opts[1]));
1387     PetscCall(PetscViewerASCIIPrintf(viewer,"        relaxation type %" PetscInt_FMT "\n",jac->as_amg_beta_opts[2]));
1388     PetscCall(PetscViewerASCIIPrintf(viewer,"        interpolation type %" PetscInt_FMT "\n",jac->as_amg_beta_opts[3]));
1389     PetscCall(PetscViewerASCIIPrintf(viewer,"        max nonzero elements in interpolation rows %" PetscInt_FMT "\n",jac->as_amg_beta_opts[4]));
1390     PetscCall(PetscViewerASCIIPrintf(viewer,"        strength threshold %g\n",(double)jac->as_amg_beta_theta));
1391   }
1392   PetscFunctionReturn(0);
1393 }
1394 
1395 static PetscErrorCode PCHYPRESetDiscreteGradient_HYPRE(PC pc, Mat G)
1396 {
1397   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1398   PetscBool      ishypre;
1399 
1400   PetscFunctionBegin;
1401   PetscCall(PetscObjectTypeCompare((PetscObject)G,MATHYPRE,&ishypre));
1402   if (ishypre) {
1403     PetscCall(PetscObjectReference((PetscObject)G));
1404     PetscCall(MatDestroy(&jac->G));
1405     jac->G = G;
1406   } else {
1407     PetscCall(MatDestroy(&jac->G));
1408     PetscCall(MatConvert(G,MATHYPRE,MAT_INITIAL_MATRIX,&jac->G));
1409   }
1410   PetscFunctionReturn(0);
1411 }
1412 
1413 /*@
1414  PCHYPRESetDiscreteGradient - Set discrete gradient matrix
1415 
1416    Collective on PC
1417 
1418    Input Parameters:
1419 +  pc - the preconditioning context
1420 -  G - the discrete gradient
1421 
1422    Level: intermediate
1423 
1424    Notes:
1425     G should have as many rows as the number of edges and as many columns as the number of vertices in the mesh
1426 
1427     Each row of G has 2 nonzeros, with column indexes being the global indexes of edge's endpoints: matrix entries are +1 and -1 depending on edge orientation
1428 
1429 .seealso: `PCHYPRESetDiscreteCurl()`
1430 @*/
1431 PetscErrorCode PCHYPRESetDiscreteGradient(PC pc, Mat G)
1432 {
1433   PetscFunctionBegin;
1434   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1435   PetscValidHeaderSpecific(G,MAT_CLASSID,2);
1436   PetscCheckSameComm(pc,1,G,2);
1437   PetscTryMethod(pc,"PCHYPRESetDiscreteGradient_C",(PC,Mat),(pc,G));
1438   PetscFunctionReturn(0);
1439 }
1440 
1441 static PetscErrorCode PCHYPRESetDiscreteCurl_HYPRE(PC pc, Mat C)
1442 {
1443   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1444   PetscBool      ishypre;
1445 
1446   PetscFunctionBegin;
1447   PetscCall(PetscObjectTypeCompare((PetscObject)C,MATHYPRE,&ishypre));
1448   if (ishypre) {
1449     PetscCall(PetscObjectReference((PetscObject)C));
1450     PetscCall(MatDestroy(&jac->C));
1451     jac->C = C;
1452   } else {
1453     PetscCall(MatDestroy(&jac->C));
1454     PetscCall(MatConvert(C,MATHYPRE,MAT_INITIAL_MATRIX,&jac->C));
1455   }
1456   PetscFunctionReturn(0);
1457 }
1458 
1459 /*@
1460  PCHYPRESetDiscreteCurl - Set discrete curl matrix
1461 
1462    Collective on PC
1463 
1464    Input Parameters:
1465 +  pc - the preconditioning context
1466 -  C - the discrete curl
1467 
1468    Level: intermediate
1469 
1470    Notes:
1471     C should have as many rows as the number of faces and as many columns as the number of edges in the mesh
1472 
1473     Each row of G has as many nonzeros as the number of edges of a face, with column indexes being the global indexes of the corresponding edge: matrix entries are +1 and -1 depending on edge orientation with respect to the face orientation
1474 
1475 .seealso: `PCHYPRESetDiscreteGradient()`
1476 @*/
1477 PetscErrorCode PCHYPRESetDiscreteCurl(PC pc, Mat C)
1478 {
1479   PetscFunctionBegin;
1480   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1481   PetscValidHeaderSpecific(C,MAT_CLASSID,2);
1482   PetscCheckSameComm(pc,1,C,2);
1483   PetscTryMethod(pc,"PCHYPRESetDiscreteCurl_C",(PC,Mat),(pc,C));
1484   PetscFunctionReturn(0);
1485 }
1486 
1487 static PetscErrorCode PCHYPRESetInterpolations_HYPRE(PC pc, PetscInt dim, Mat RT_PiFull, Mat RT_Pi[], Mat ND_PiFull, Mat ND_Pi[])
1488 {
1489   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1490   PetscBool      ishypre;
1491   PetscInt       i;
1492   PetscFunctionBegin;
1493 
1494   PetscCall(MatDestroy(&jac->RT_PiFull));
1495   PetscCall(MatDestroy(&jac->ND_PiFull));
1496   for (i=0;i<3;++i) {
1497     PetscCall(MatDestroy(&jac->RT_Pi[i]));
1498     PetscCall(MatDestroy(&jac->ND_Pi[i]));
1499   }
1500 
1501   jac->dim = dim;
1502   if (RT_PiFull) {
1503     PetscCall(PetscObjectTypeCompare((PetscObject)RT_PiFull,MATHYPRE,&ishypre));
1504     if (ishypre) {
1505       PetscCall(PetscObjectReference((PetscObject)RT_PiFull));
1506       jac->RT_PiFull = RT_PiFull;
1507     } else {
1508       PetscCall(MatConvert(RT_PiFull,MATHYPRE,MAT_INITIAL_MATRIX,&jac->RT_PiFull));
1509     }
1510   }
1511   if (RT_Pi) {
1512     for (i=0;i<dim;++i) {
1513       if (RT_Pi[i]) {
1514         PetscCall(PetscObjectTypeCompare((PetscObject)RT_Pi[i],MATHYPRE,&ishypre));
1515         if (ishypre) {
1516           PetscCall(PetscObjectReference((PetscObject)RT_Pi[i]));
1517           jac->RT_Pi[i] = RT_Pi[i];
1518         } else {
1519           PetscCall(MatConvert(RT_Pi[i],MATHYPRE,MAT_INITIAL_MATRIX,&jac->RT_Pi[i]));
1520         }
1521       }
1522     }
1523   }
1524   if (ND_PiFull) {
1525     PetscCall(PetscObjectTypeCompare((PetscObject)ND_PiFull,MATHYPRE,&ishypre));
1526     if (ishypre) {
1527       PetscCall(PetscObjectReference((PetscObject)ND_PiFull));
1528       jac->ND_PiFull = ND_PiFull;
1529     } else {
1530       PetscCall(MatConvert(ND_PiFull,MATHYPRE,MAT_INITIAL_MATRIX,&jac->ND_PiFull));
1531     }
1532   }
1533   if (ND_Pi) {
1534     for (i=0;i<dim;++i) {
1535       if (ND_Pi[i]) {
1536         PetscCall(PetscObjectTypeCompare((PetscObject)ND_Pi[i],MATHYPRE,&ishypre));
1537         if (ishypre) {
1538           PetscCall(PetscObjectReference((PetscObject)ND_Pi[i]));
1539           jac->ND_Pi[i] = ND_Pi[i];
1540         } else {
1541           PetscCall(MatConvert(ND_Pi[i],MATHYPRE,MAT_INITIAL_MATRIX,&jac->ND_Pi[i]));
1542         }
1543       }
1544     }
1545   }
1546 
1547   PetscFunctionReturn(0);
1548 }
1549 
1550 /*@
1551  PCHYPRESetInterpolations - Set interpolation matrices for AMS/ADS preconditioner
1552 
1553    Collective on PC
1554 
1555    Input Parameters:
1556 +  pc - the preconditioning context
1557 -  dim - the dimension of the problem, only used in AMS
1558 -  RT_PiFull - Raviart-Thomas interpolation matrix
1559 -  RT_Pi - x/y/z component of Raviart-Thomas interpolation matrix
1560 -  ND_PiFull - Nedelec interpolation matrix
1561 -  ND_Pi - x/y/z component of Nedelec interpolation matrix
1562 
1563    Notes:
1564     For AMS, only Nedelec interpolation matrices are needed, the Raviart-Thomas interpolation matrices can be set to NULL.
1565 
1566     For ADS, both type of interpolation matrices are needed.
1567 
1568    Level: intermediate
1569 
1570 @*/
1571 PetscErrorCode PCHYPRESetInterpolations(PC pc, PetscInt dim, Mat RT_PiFull, Mat RT_Pi[], Mat ND_PiFull, Mat ND_Pi[])
1572 {
1573   PetscInt       i;
1574 
1575   PetscFunctionBegin;
1576   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1577   if (RT_PiFull) {
1578     PetscValidHeaderSpecific(RT_PiFull,MAT_CLASSID,3);
1579     PetscCheckSameComm(pc,1,RT_PiFull,3);
1580   }
1581   if (RT_Pi) {
1582     PetscValidPointer(RT_Pi,4);
1583     for (i=0;i<dim;++i) {
1584       if (RT_Pi[i]) {
1585         PetscValidHeaderSpecific(RT_Pi[i],MAT_CLASSID,4);
1586         PetscCheckSameComm(pc,1,RT_Pi[i],4);
1587       }
1588     }
1589   }
1590   if (ND_PiFull) {
1591     PetscValidHeaderSpecific(ND_PiFull,MAT_CLASSID,5);
1592     PetscCheckSameComm(pc,1,ND_PiFull,5);
1593   }
1594   if (ND_Pi) {
1595     PetscValidPointer(ND_Pi,6);
1596     for (i=0;i<dim;++i) {
1597       if (ND_Pi[i]) {
1598         PetscValidHeaderSpecific(ND_Pi[i],MAT_CLASSID,6);
1599         PetscCheckSameComm(pc,1,ND_Pi[i],6);
1600       }
1601     }
1602   }
1603   PetscTryMethod(pc,"PCHYPRESetInterpolations_C",(PC,PetscInt,Mat,Mat[],Mat,Mat[]),(pc,dim,RT_PiFull,RT_Pi,ND_PiFull,ND_Pi));
1604   PetscFunctionReturn(0);
1605 }
1606 
1607 static PetscErrorCode PCHYPRESetPoissonMatrix_HYPRE(PC pc, Mat A, PetscBool isalpha)
1608 {
1609   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1610   PetscBool      ishypre;
1611 
1612   PetscFunctionBegin;
1613   PetscCall(PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre));
1614   if (ishypre) {
1615     if (isalpha) {
1616       PetscCall(PetscObjectReference((PetscObject)A));
1617       PetscCall(MatDestroy(&jac->alpha_Poisson));
1618       jac->alpha_Poisson = A;
1619     } else {
1620       if (A) {
1621         PetscCall(PetscObjectReference((PetscObject)A));
1622       } else {
1623         jac->ams_beta_is_zero = PETSC_TRUE;
1624       }
1625       PetscCall(MatDestroy(&jac->beta_Poisson));
1626       jac->beta_Poisson = A;
1627     }
1628   } else {
1629     if (isalpha) {
1630       PetscCall(MatDestroy(&jac->alpha_Poisson));
1631       PetscCall(MatConvert(A,MATHYPRE,MAT_INITIAL_MATRIX,&jac->alpha_Poisson));
1632     } else {
1633       if (A) {
1634         PetscCall(MatDestroy(&jac->beta_Poisson));
1635         PetscCall(MatConvert(A,MATHYPRE,MAT_INITIAL_MATRIX,&jac->beta_Poisson));
1636       } else {
1637         PetscCall(MatDestroy(&jac->beta_Poisson));
1638         jac->ams_beta_is_zero = PETSC_TRUE;
1639       }
1640     }
1641   }
1642   PetscFunctionReturn(0);
1643 }
1644 
1645 /*@
1646  PCHYPRESetAlphaPoissonMatrix - Set vector Poisson matrix
1647 
1648    Collective on PC
1649 
1650    Input Parameters:
1651 +  pc - the preconditioning context
1652 -  A - the matrix
1653 
1654    Level: intermediate
1655 
1656    Notes:
1657     A should be obtained by discretizing the vector valued Poisson problem with linear finite elements
1658 
1659 .seealso: `PCHYPRESetDiscreteGradient()`, `PCHYPRESetDiscreteCurl()`, `PCHYPRESetBetaPoissonMatrix()`
1660 @*/
1661 PetscErrorCode PCHYPRESetAlphaPoissonMatrix(PC pc, Mat A)
1662 {
1663   PetscFunctionBegin;
1664   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1665   PetscValidHeaderSpecific(A,MAT_CLASSID,2);
1666   PetscCheckSameComm(pc,1,A,2);
1667   PetscTryMethod(pc,"PCHYPRESetPoissonMatrix_C",(PC,Mat,PetscBool),(pc,A,PETSC_TRUE));
1668   PetscFunctionReturn(0);
1669 }
1670 
1671 /*@
1672  PCHYPRESetBetaPoissonMatrix - Set Poisson matrix
1673 
1674    Collective on PC
1675 
1676    Input Parameters:
1677 +  pc - the preconditioning context
1678 -  A - the matrix
1679 
1680    Level: intermediate
1681 
1682    Notes:
1683     A should be obtained by discretizing the Poisson problem with linear finite elements.
1684           Following HYPRE convention, the scalar Poisson solver of AMS can be turned off by passing NULL.
1685 
1686 .seealso: `PCHYPRESetDiscreteGradient()`, `PCHYPRESetDiscreteCurl()`, `PCHYPRESetAlphaPoissonMatrix()`
1687 @*/
1688 PetscErrorCode PCHYPRESetBetaPoissonMatrix(PC pc, Mat A)
1689 {
1690   PetscFunctionBegin;
1691   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1692   if (A) {
1693     PetscValidHeaderSpecific(A,MAT_CLASSID,2);
1694     PetscCheckSameComm(pc,1,A,2);
1695   }
1696   PetscTryMethod(pc,"PCHYPRESetPoissonMatrix_C",(PC,Mat,PetscBool),(pc,A,PETSC_FALSE));
1697   PetscFunctionReturn(0);
1698 }
1699 
1700 static PetscErrorCode PCHYPRESetEdgeConstantVectors_HYPRE(PC pc,Vec ozz, Vec zoz, Vec zzo)
1701 {
1702   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
1703 
1704   PetscFunctionBegin;
1705   /* throw away any vector if already set */
1706   PetscCall(VecHYPRE_IJVectorDestroy(&jac->constants[0]));
1707   PetscCall(VecHYPRE_IJVectorDestroy(&jac->constants[1]));
1708   PetscCall(VecHYPRE_IJVectorDestroy(&jac->constants[2]));
1709   PetscCall(VecHYPRE_IJVectorCreate(ozz->map,&jac->constants[0]));
1710   PetscCall(VecHYPRE_IJVectorCopy(ozz,jac->constants[0]));
1711   PetscCall(VecHYPRE_IJVectorCreate(zoz->map,&jac->constants[1]));
1712   PetscCall(VecHYPRE_IJVectorCopy(zoz,jac->constants[1]));
1713   jac->dim = 2;
1714   if (zzo) {
1715     PetscCall(VecHYPRE_IJVectorCreate(zzo->map,&jac->constants[2]));
1716     PetscCall(VecHYPRE_IJVectorCopy(zzo,jac->constants[2]));
1717     jac->dim++;
1718   }
1719   PetscFunctionReturn(0);
1720 }
1721 
1722 /*@
1723  PCHYPRESetEdgeConstantVectors - Set the representation of the constant vector fields in the edge element basis
1724 
1725    Collective on PC
1726 
1727    Input Parameters:
1728 +  pc - the preconditioning context
1729 -  ozz - vector representing (1,0,0) (or (1,0) in 2D)
1730 -  zoz - vector representing (0,1,0) (or (0,1) in 2D)
1731 -  zzo - vector representing (0,0,1) (use NULL in 2D)
1732 
1733    Level: intermediate
1734 
1735 @*/
1736 PetscErrorCode PCHYPRESetEdgeConstantVectors(PC pc, Vec ozz, Vec zoz, Vec zzo)
1737 {
1738   PetscFunctionBegin;
1739   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1740   PetscValidHeaderSpecific(ozz,VEC_CLASSID,2);
1741   PetscValidHeaderSpecific(zoz,VEC_CLASSID,3);
1742   if (zzo) PetscValidHeaderSpecific(zzo,VEC_CLASSID,4);
1743   PetscCheckSameComm(pc,1,ozz,2);
1744   PetscCheckSameComm(pc,1,zoz,3);
1745   if (zzo) PetscCheckSameComm(pc,1,zzo,4);
1746   PetscTryMethod(pc,"PCHYPRESetEdgeConstantVectors_C",(PC,Vec,Vec,Vec),(pc,ozz,zoz,zzo));
1747   PetscFunctionReturn(0);
1748 }
1749 
1750 static PetscErrorCode PCHYPREAMSSetInteriorNodes_HYPRE(PC pc,Vec interior)
1751 {
1752   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
1753 
1754   PetscFunctionBegin;
1755   PetscCall(VecHYPRE_IJVectorDestroy(&jac->interior));
1756   PetscCall(VecHYPRE_IJVectorCreate(interior->map,&jac->interior));
1757   PetscCall(VecHYPRE_IJVectorCopy(interior,jac->interior));
1758   jac->ams_beta_is_zero_part = PETSC_TRUE;
1759   PetscFunctionReturn(0);
1760 }
1761 
1762 /*@
1763  PCHYPREAMSSetInteriorNodes - Set the list of interior nodes to a zero-conductivity region.
1764 
1765    Collective on PC
1766 
1767    Input Parameters:
1768 +  pc - the preconditioning context
1769 -  interior - vector. node is interior if its entry in the array is 1.0.
1770 
1771    Level: intermediate
1772 
1773    Note:
1774    This calls HYPRE_AMSSetInteriorNodes()
1775 .seealso:
1776 @*/
1777 PetscErrorCode PCHYPREAMSSetInteriorNodes(PC pc, Vec interior)
1778 {
1779   PetscFunctionBegin;
1780   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1781   PetscValidHeaderSpecific(interior,VEC_CLASSID,2);
1782   PetscCheckSameComm(pc,1,interior,2);
1783   PetscTryMethod(pc,"PCHYPREAMSSetInteriorNodes_C",(PC,Vec),(pc,interior));
1784   PetscFunctionReturn(0);
1785 }
1786 
1787 static PetscErrorCode PCSetCoordinates_HYPRE(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords)
1788 {
1789   PC_HYPRE        *jac = (PC_HYPRE*)pc->data;
1790   Vec             tv;
1791   PetscInt        i;
1792 
1793   PetscFunctionBegin;
1794   /* throw away any coordinate vector if already set */
1795   PetscCall(VecHYPRE_IJVectorDestroy(&jac->coords[0]));
1796   PetscCall(VecHYPRE_IJVectorDestroy(&jac->coords[1]));
1797   PetscCall(VecHYPRE_IJVectorDestroy(&jac->coords[2]));
1798   jac->dim = dim;
1799 
1800   /* compute IJ vector for coordinates */
1801   PetscCall(VecCreate(PetscObjectComm((PetscObject)pc),&tv));
1802   PetscCall(VecSetType(tv,VECSTANDARD));
1803   PetscCall(VecSetSizes(tv,nloc,PETSC_DECIDE));
1804   for (i=0;i<dim;i++) {
1805     PetscScalar *array;
1806     PetscInt    j;
1807 
1808     PetscCall(VecHYPRE_IJVectorCreate(tv->map,&jac->coords[i]));
1809     PetscCall(VecGetArrayWrite(tv,&array));
1810     for (j=0;j<nloc;j++) array[j] = coords[j*dim+i];
1811     PetscCall(VecRestoreArrayWrite(tv,&array));
1812     PetscCall(VecHYPRE_IJVectorCopy(tv,jac->coords[i]));
1813   }
1814   PetscCall(VecDestroy(&tv));
1815   PetscFunctionReturn(0);
1816 }
1817 
1818 /* ---------------------------------------------------------------------------------*/
1819 
1820 static PetscErrorCode  PCHYPREGetType_HYPRE(PC pc,const char *name[])
1821 {
1822   PC_HYPRE *jac = (PC_HYPRE*)pc->data;
1823 
1824   PetscFunctionBegin;
1825   *name = jac->hypre_type;
1826   PetscFunctionReturn(0);
1827 }
1828 
1829 static PetscErrorCode  PCHYPRESetType_HYPRE(PC pc,const char name[])
1830 {
1831   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1832   PetscBool      flag;
1833 
1834   PetscFunctionBegin;
1835   if (jac->hypre_type) {
1836     PetscCall(PetscStrcmp(jac->hypre_type,name,&flag));
1837     PetscCheck(flag,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
1838     PetscFunctionReturn(0);
1839   } else {
1840     PetscCall(PetscStrallocpy(name, &jac->hypre_type));
1841   }
1842 
1843   jac->maxiter         = PETSC_DEFAULT;
1844   jac->tol             = PETSC_DEFAULT;
1845   jac->printstatistics = PetscLogPrintInfo;
1846 
1847   PetscCall(PetscStrcmp("pilut",jac->hypre_type,&flag));
1848   if (flag) {
1849     PetscCall(PetscCommGetComm(PetscObjectComm((PetscObject)pc),&jac->comm_hypre));
1850     PetscCallExternal(HYPRE_ParCSRPilutCreate,jac->comm_hypre,&jac->hsolver);
1851     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
1852     pc->ops->view           = PCView_HYPRE_Pilut;
1853     jac->destroy            = HYPRE_ParCSRPilutDestroy;
1854     jac->setup              = HYPRE_ParCSRPilutSetup;
1855     jac->solve              = HYPRE_ParCSRPilutSolve;
1856     jac->factorrowsize      = PETSC_DEFAULT;
1857     PetscFunctionReturn(0);
1858   }
1859   PetscCall(PetscStrcmp("euclid",jac->hypre_type,&flag));
1860   if (flag) {
1861 #if defined(PETSC_USE_64BIT_INDICES)
1862     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Hypre Euclid does not support 64 bit indices");
1863 #endif
1864     PetscCall(PetscCommGetComm(PetscObjectComm((PetscObject)pc),&jac->comm_hypre));
1865     PetscCallExternal(HYPRE_EuclidCreate,jac->comm_hypre,&jac->hsolver);
1866     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
1867     pc->ops->view           = PCView_HYPRE_Euclid;
1868     jac->destroy            = HYPRE_EuclidDestroy;
1869     jac->setup              = HYPRE_EuclidSetup;
1870     jac->solve              = HYPRE_EuclidSolve;
1871     jac->factorrowsize      = PETSC_DEFAULT;
1872     jac->eu_level           = PETSC_DEFAULT; /* default */
1873     PetscFunctionReturn(0);
1874   }
1875   PetscCall(PetscStrcmp("parasails",jac->hypre_type,&flag));
1876   if (flag) {
1877     PetscCall(PetscCommGetComm(PetscObjectComm((PetscObject)pc),&jac->comm_hypre));
1878     PetscCallExternal(HYPRE_ParaSailsCreate,jac->comm_hypre,&jac->hsolver);
1879     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
1880     pc->ops->view           = PCView_HYPRE_ParaSails;
1881     jac->destroy            = HYPRE_ParaSailsDestroy;
1882     jac->setup              = HYPRE_ParaSailsSetup;
1883     jac->solve              = HYPRE_ParaSailsSolve;
1884     /* initialize */
1885     jac->nlevels   = 1;
1886     jac->threshold = .1;
1887     jac->filter    = .1;
1888     jac->loadbal   = 0;
1889     if (PetscLogPrintInfo) jac->logging = (int) PETSC_TRUE;
1890     else jac->logging = (int) PETSC_FALSE;
1891 
1892     jac->ruse = (int) PETSC_FALSE;
1893     jac->symt = 0;
1894     PetscCallExternal(HYPRE_ParaSailsSetParams,jac->hsolver,jac->threshold,jac->nlevels);
1895     PetscCallExternal(HYPRE_ParaSailsSetFilter,jac->hsolver,jac->filter);
1896     PetscCallExternal(HYPRE_ParaSailsSetLoadbal,jac->hsolver,jac->loadbal);
1897     PetscCallExternal(HYPRE_ParaSailsSetLogging,jac->hsolver,jac->logging);
1898     PetscCallExternal(HYPRE_ParaSailsSetReuse,jac->hsolver,jac->ruse);
1899     PetscCallExternal(HYPRE_ParaSailsSetSym,jac->hsolver,jac->symt);
1900     PetscFunctionReturn(0);
1901   }
1902   PetscCall(PetscStrcmp("boomeramg",jac->hypre_type,&flag));
1903   if (flag) {
1904     PetscCallExternal(HYPRE_BoomerAMGCreate,&jac->hsolver);
1905     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
1906     pc->ops->view            = PCView_HYPRE_BoomerAMG;
1907     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
1908     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
1909     PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_BoomerAMG));
1910     PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_BoomerAMG));
1911     jac->destroy             = HYPRE_BoomerAMGDestroy;
1912     jac->setup               = HYPRE_BoomerAMGSetup;
1913     jac->solve               = HYPRE_BoomerAMGSolve;
1914     jac->applyrichardson     = PETSC_FALSE;
1915     /* these defaults match the hypre defaults */
1916     jac->cycletype        = 1;
1917     jac->maxlevels        = 25;
1918     jac->maxiter          = 1;
1919     jac->tol              = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */
1920     jac->truncfactor      = 0.0;
1921     jac->strongthreshold  = .25;
1922     jac->maxrowsum        = .9;
1923     jac->coarsentype      = 6;
1924     jac->measuretype      = 0;
1925     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = 1;
1926     jac->smoothtype       = -1; /* Not set by default */
1927     jac->smoothnumlevels  = 25;
1928     jac->eu_level         = 0;
1929     jac->eu_droptolerance = 0;
1930     jac->eu_bj            = 0;
1931     jac->relaxtype[0]     = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a PC - most likely with CG */
1932     jac->relaxtype[2]     = 9; /*G.E. */
1933     jac->relaxweight      = 1.0;
1934     jac->outerrelaxweight = 1.0;
1935     jac->relaxorder       = 1;
1936     jac->interptype       = 0;
1937     jac->Rtype            = 0;
1938     jac->Rstrongthreshold = 0.25;
1939     jac->Rfilterthreshold = 0.0;
1940     jac->Adroptype        = -1;
1941     jac->Adroptol         = 0.0;
1942     jac->agg_nl           = 0;
1943     jac->agg_interptype   = 4;
1944     jac->pmax             = 0;
1945     jac->truncfactor      = 0.0;
1946     jac->agg_num_paths    = 1;
1947     jac->maxc             = 9;
1948     jac->minc             = 1;
1949     jac->nodal_coarsening      = 0;
1950     jac->nodal_coarsening_diag = 0;
1951     jac->vec_interp_variant    = 0;
1952     jac->vec_interp_qmax       = 0;
1953     jac->vec_interp_smooth     = PETSC_FALSE;
1954     jac->interp_refine         = 0;
1955     jac->nodal_relax           = PETSC_FALSE;
1956     jac->nodal_relax_levels    = 1;
1957     jac->rap2                  = 0;
1958 
1959     /* GPU defaults
1960          from https://hypre.readthedocs.io/en/latest/solvers-boomeramg.html#gpu-supported-options
1961          and /src/parcsr_ls/par_amg.c */
1962 #if defined(PETSC_HAVE_HYPRE_DEVICE)
1963     jac->keeptranspose         = PETSC_TRUE;
1964     jac->mod_rap2              = 1;
1965     jac->coarsentype           = 8;
1966     jac->relaxorder            = 0;
1967     jac->interptype            = 6;
1968     jac->relaxtype[0]          = 18;
1969     jac->relaxtype[1]          = 18;
1970     jac->agg_interptype        = 7;
1971 #else
1972     jac->keeptranspose         = PETSC_FALSE;
1973     jac->mod_rap2              = 0;
1974 #endif
1975     PetscCallExternal(HYPRE_BoomerAMGSetCycleType,jac->hsolver,jac->cycletype);
1976     PetscCallExternal(HYPRE_BoomerAMGSetMaxLevels,jac->hsolver,jac->maxlevels);
1977     PetscCallExternal(HYPRE_BoomerAMGSetMaxIter,jac->hsolver,jac->maxiter);
1978     PetscCallExternal(HYPRE_BoomerAMGSetTol,jac->hsolver,jac->tol);
1979     PetscCallExternal(HYPRE_BoomerAMGSetTruncFactor,jac->hsolver,jac->truncfactor);
1980     PetscCallExternal(HYPRE_BoomerAMGSetStrongThreshold,jac->hsolver,jac->strongthreshold);
1981     PetscCallExternal(HYPRE_BoomerAMGSetMaxRowSum,jac->hsolver,jac->maxrowsum);
1982     PetscCallExternal(HYPRE_BoomerAMGSetCoarsenType,jac->hsolver,jac->coarsentype);
1983     PetscCallExternal(HYPRE_BoomerAMGSetMeasureType,jac->hsolver,jac->measuretype);
1984     PetscCallExternal(HYPRE_BoomerAMGSetRelaxOrder,jac->hsolver, jac->relaxorder);
1985     PetscCallExternal(HYPRE_BoomerAMGSetInterpType,jac->hsolver,jac->interptype);
1986     PetscCallExternal(HYPRE_BoomerAMGSetAggNumLevels,jac->hsolver,jac->agg_nl);
1987     PetscCallExternal(HYPRE_BoomerAMGSetAggInterpType,jac->hsolver,jac->agg_interptype);
1988     PetscCallExternal(HYPRE_BoomerAMGSetPMaxElmts,jac->hsolver,jac->pmax);
1989     PetscCallExternal(HYPRE_BoomerAMGSetNumPaths,jac->hsolver,jac->agg_num_paths);
1990     PetscCallExternal(HYPRE_BoomerAMGSetRelaxType,jac->hsolver, jac->relaxtype[0]);  /* defaults coarse to 9 */
1991     PetscCallExternal(HYPRE_BoomerAMGSetNumSweeps,jac->hsolver, jac->gridsweeps[0]); /* defaults coarse to 1 */
1992     PetscCallExternal(HYPRE_BoomerAMGSetMaxCoarseSize,jac->hsolver, jac->maxc);
1993     PetscCallExternal(HYPRE_BoomerAMGSetMinCoarseSize,jac->hsolver, jac->minc);
1994     /* GPU */
1995 #if PETSC_PKG_HYPRE_VERSION_GE(2,18,0)
1996     PetscCallExternal(HYPRE_BoomerAMGSetKeepTranspose,jac->hsolver,jac->keeptranspose ? 1 : 0);
1997     PetscCallExternal(HYPRE_BoomerAMGSetRAP2,jac->hsolver, jac->rap2);
1998     PetscCallExternal(HYPRE_BoomerAMGSetModuleRAP2,jac->hsolver, jac->mod_rap2);
1999 #endif
2000 
2001     /* AIR */
2002 #if PETSC_PKG_HYPRE_VERSION_GE(2,18,0)
2003     PetscCallExternal(HYPRE_BoomerAMGSetRestriction,jac->hsolver,jac->Rtype);
2004     PetscCallExternal(HYPRE_BoomerAMGSetStrongThresholdR,jac->hsolver,jac->Rstrongthreshold);
2005     PetscCallExternal(HYPRE_BoomerAMGSetFilterThresholdR,jac->hsolver,jac->Rfilterthreshold);
2006     PetscCallExternal(HYPRE_BoomerAMGSetADropTol,jac->hsolver,jac->Adroptol);
2007     PetscCallExternal(HYPRE_BoomerAMGSetADropType,jac->hsolver,jac->Adroptype);
2008 #endif
2009     PetscFunctionReturn(0);
2010   }
2011   PetscCall(PetscStrcmp("ams",jac->hypre_type,&flag));
2012   if (flag) {
2013     PetscCall(HYPRE_AMSCreate(&jac->hsolver));
2014     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_AMS;
2015     pc->ops->view            = PCView_HYPRE_AMS;
2016     jac->destroy             = HYPRE_AMSDestroy;
2017     jac->setup               = HYPRE_AMSSetup;
2018     jac->solve               = HYPRE_AMSSolve;
2019     jac->coords[0]           = NULL;
2020     jac->coords[1]           = NULL;
2021     jac->coords[2]           = NULL;
2022     jac->interior            = NULL;
2023     /* solver parameters: these are borrowed from mfem package, and they are not the default values from HYPRE AMS */
2024     jac->as_print           = 0;
2025     jac->as_max_iter        = 1; /* used as a preconditioner */
2026     jac->as_tol             = 0.; /* used as a preconditioner */
2027     jac->ams_cycle_type     = 13;
2028     /* Smoothing options */
2029     jac->as_relax_type      = 2;
2030     jac->as_relax_times     = 1;
2031     jac->as_relax_weight    = 1.0;
2032     jac->as_omega           = 1.0;
2033     /* Vector valued Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
2034     jac->as_amg_alpha_opts[0] = 10;
2035     jac->as_amg_alpha_opts[1] = 1;
2036     jac->as_amg_alpha_opts[2] = 6;
2037     jac->as_amg_alpha_opts[3] = 6;
2038     jac->as_amg_alpha_opts[4] = 4;
2039     jac->as_amg_alpha_theta   = 0.25;
2040     /* Scalar Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
2041     jac->as_amg_beta_opts[0] = 10;
2042     jac->as_amg_beta_opts[1] = 1;
2043     jac->as_amg_beta_opts[2] = 6;
2044     jac->as_amg_beta_opts[3] = 6;
2045     jac->as_amg_beta_opts[4] = 4;
2046     jac->as_amg_beta_theta   = 0.25;
2047     PetscCallExternal(HYPRE_AMSSetPrintLevel,jac->hsolver,jac->as_print);
2048     PetscCallExternal(HYPRE_AMSSetMaxIter,jac->hsolver,jac->as_max_iter);
2049     PetscCallExternal(HYPRE_AMSSetCycleType,jac->hsolver,jac->ams_cycle_type);
2050     PetscCallExternal(HYPRE_AMSSetTol,jac->hsolver,jac->as_tol);
2051     PetscCallExternal(HYPRE_AMSSetSmoothingOptions,jac->hsolver,jac->as_relax_type,
2052                                                                       jac->as_relax_times,
2053                                                                       jac->as_relax_weight,
2054                                                                       jac->as_omega);
2055     PetscCallExternal(HYPRE_AMSSetAlphaAMGOptions,jac->hsolver,jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
2056                                                                      jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
2057                                                                      jac->as_amg_alpha_opts[2],       /* AMG relax_type */
2058                                                                      jac->as_amg_alpha_theta,
2059                                                                      jac->as_amg_alpha_opts[3],       /* AMG interp_type */
2060                                                                      jac->as_amg_alpha_opts[4]);     /* AMG Pmax */
2061     PetscCallExternal(HYPRE_AMSSetBetaAMGOptions,jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
2062                                                                     jac->as_amg_beta_opts[1],       /* AMG agg_levels */
2063                                                                     jac->as_amg_beta_opts[2],       /* AMG relax_type */
2064                                                                     jac->as_amg_beta_theta,
2065                                                                     jac->as_amg_beta_opts[3],       /* AMG interp_type */
2066                                                                     jac->as_amg_beta_opts[4]);     /* AMG Pmax */
2067     /* Zero conductivity */
2068     jac->ams_beta_is_zero      = PETSC_FALSE;
2069     jac->ams_beta_is_zero_part = PETSC_FALSE;
2070     PetscFunctionReturn(0);
2071   }
2072   PetscCall(PetscStrcmp("ads",jac->hypre_type,&flag));
2073   if (flag) {
2074     PetscCall(HYPRE_ADSCreate(&jac->hsolver));
2075     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_ADS;
2076     pc->ops->view            = PCView_HYPRE_ADS;
2077     jac->destroy             = HYPRE_ADSDestroy;
2078     jac->setup               = HYPRE_ADSSetup;
2079     jac->solve               = HYPRE_ADSSolve;
2080     jac->coords[0]           = NULL;
2081     jac->coords[1]           = NULL;
2082     jac->coords[2]           = NULL;
2083     /* solver parameters: these are borrowed from mfem package, and they are not the default values from HYPRE ADS */
2084     jac->as_print           = 0;
2085     jac->as_max_iter        = 1; /* used as a preconditioner */
2086     jac->as_tol             = 0.; /* used as a preconditioner */
2087     jac->ads_cycle_type     = 13;
2088     /* Smoothing options */
2089     jac->as_relax_type      = 2;
2090     jac->as_relax_times     = 1;
2091     jac->as_relax_weight    = 1.0;
2092     jac->as_omega           = 1.0;
2093     /* AMS solver parameters: cycle_type, coarsen type, agg_levels, relax_type, interp_type, Pmax */
2094     jac->ams_cycle_type       = 14;
2095     jac->as_amg_alpha_opts[0] = 10;
2096     jac->as_amg_alpha_opts[1] = 1;
2097     jac->as_amg_alpha_opts[2] = 6;
2098     jac->as_amg_alpha_opts[3] = 6;
2099     jac->as_amg_alpha_opts[4] = 4;
2100     jac->as_amg_alpha_theta   = 0.25;
2101     /* Vector Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
2102     jac->as_amg_beta_opts[0] = 10;
2103     jac->as_amg_beta_opts[1] = 1;
2104     jac->as_amg_beta_opts[2] = 6;
2105     jac->as_amg_beta_opts[3] = 6;
2106     jac->as_amg_beta_opts[4] = 4;
2107     jac->as_amg_beta_theta   = 0.25;
2108     PetscCallExternal(HYPRE_ADSSetPrintLevel,jac->hsolver,jac->as_print);
2109     PetscCallExternal(HYPRE_ADSSetMaxIter,jac->hsolver,jac->as_max_iter);
2110     PetscCallExternal(HYPRE_ADSSetCycleType,jac->hsolver,jac->ams_cycle_type);
2111     PetscCallExternal(HYPRE_ADSSetTol,jac->hsolver,jac->as_tol);
2112     PetscCallExternal(HYPRE_ADSSetSmoothingOptions,jac->hsolver,jac->as_relax_type,
2113                                                                       jac->as_relax_times,
2114                                                                       jac->as_relax_weight,
2115                                                                       jac->as_omega);
2116     PetscCallExternal(HYPRE_ADSSetAMSOptions,jac->hsolver,jac->ams_cycle_type,             /* AMG coarsen type */
2117                                                                 jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
2118                                                                 jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
2119                                                                 jac->as_amg_alpha_opts[2],       /* AMG relax_type */
2120                                                                 jac->as_amg_alpha_theta,
2121                                                                 jac->as_amg_alpha_opts[3],       /* AMG interp_type */
2122                                                                 jac->as_amg_alpha_opts[4]);     /* AMG Pmax */
2123     PetscCallExternal(HYPRE_ADSSetAMGOptions,jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
2124                                                                 jac->as_amg_beta_opts[1],       /* AMG agg_levels */
2125                                                                 jac->as_amg_beta_opts[2],       /* AMG relax_type */
2126                                                                 jac->as_amg_beta_theta,
2127                                                                 jac->as_amg_beta_opts[3],       /* AMG interp_type */
2128                                                                 jac->as_amg_beta_opts[4]);     /* AMG Pmax */
2129     PetscFunctionReturn(0);
2130   }
2131   PetscCall(PetscFree(jac->hypre_type));
2132 
2133   jac->hypre_type = NULL;
2134   SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are euclid, pilut, parasails, boomeramg, ams",name);
2135 }
2136 
2137 /*
2138     It only gets here if the HYPRE type has not been set before the call to
2139    ...SetFromOptions() which actually is most of the time
2140 */
2141 PetscErrorCode PCSetFromOptions_HYPRE(PetscOptionItems *PetscOptionsObject,PC pc)
2142 {
2143   PetscInt       indx;
2144   const char     *type[] = {"euclid","pilut","parasails","boomeramg","ams","ads"};
2145   PetscBool      flg;
2146 
2147   PetscFunctionBegin;
2148   PetscOptionsHeadBegin(PetscOptionsObject,"HYPRE preconditioner options");
2149   PetscCall(PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,PETSC_STATIC_ARRAY_LENGTH(type),"boomeramg",&indx,&flg));
2150   if (flg) {
2151     PetscCall(PCHYPRESetType_HYPRE(pc,type[indx]));
2152   } else {
2153     PetscCall(PCHYPRESetType_HYPRE(pc,"boomeramg"));
2154   }
2155   if (pc->ops->setfromoptions) PetscCall(pc->ops->setfromoptions(PetscOptionsObject,pc));
2156   PetscOptionsHeadEnd();
2157   PetscFunctionReturn(0);
2158 }
2159 
2160 /*@C
2161      PCHYPRESetType - Sets which hypre preconditioner you wish to use
2162 
2163    Input Parameters:
2164 +     pc - the preconditioner context
2165 -     name - either  euclid, pilut, parasails, boomeramg, ams, ads
2166 
2167    Options Database Keys:
2168    -pc_hypre_type - One of euclid, pilut, parasails, boomeramg, ams, ads
2169 
2170    Level: intermediate
2171 
2172 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
2173           `PCHYPRE`
2174 
2175 @*/
2176 PetscErrorCode  PCHYPRESetType(PC pc,const char name[])
2177 {
2178   PetscFunctionBegin;
2179   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2180   PetscValidCharPointer(name,2);
2181   PetscTryMethod(pc,"PCHYPRESetType_C",(PC,const char[]),(pc,name));
2182   PetscFunctionReturn(0);
2183 }
2184 
2185 /*@C
2186      PCHYPREGetType - Gets which hypre preconditioner you are using
2187 
2188    Input Parameter:
2189 .     pc - the preconditioner context
2190 
2191    Output Parameter:
2192 .     name - either  euclid, pilut, parasails, boomeramg, ams, ads
2193 
2194    Level: intermediate
2195 
2196 .seealso: `PCCreate()`, `PCHYPRESetType()`, `PCType`, `PC`,
2197           `PCHYPRE`
2198 
2199 @*/
2200 PetscErrorCode  PCHYPREGetType(PC pc,const char *name[])
2201 {
2202   PetscFunctionBegin;
2203   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2204   PetscValidPointer(name,2);
2205   PetscTryMethod(pc,"PCHYPREGetType_C",(PC,const char*[]),(pc,name));
2206   PetscFunctionReturn(0);
2207 }
2208 
2209 /*@C
2210    PCMGGalerkinSetMatProductAlgorithm - Set type of SpGEMM for hypre to use
2211 
2212    Logically Collective on PC
2213 
2214    Input Parameters:
2215 +  pc - the hypre context
2216 -  type - one of 'cusparse', 'hypre'
2217 
2218    Options Database Key:
2219 .  -pc_mg_galerkin_mat_product_algorithm <cusparse,hypre> - Type of SpGEMM to use in hypre
2220 
2221    Level: intermediate
2222 
2223 .seealso: `PCMGGalerkinGetMatProductAlgorithm()`
2224 
2225 @*/
2226 PetscErrorCode PCMGGalerkinSetMatProductAlgorithm(PC pc,const char name[])
2227 {
2228   PetscFunctionBegin;
2229   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2230   PetscTryMethod(pc,"PCMGGalerkinSetMatProductAlgorithm_C",(PC,const char[]),(pc,name));
2231   PetscFunctionReturn(0);
2232 }
2233 
2234 /*@C
2235    PCMGGalerkinGetMatProductAlgorithm - Get type of SpGEMM for hypre
2236 
2237    Not Collective
2238 
2239    Input Parameter:
2240 .  pc - the multigrid context
2241 
2242    Output Parameter:
2243 .  name - one of 'cusparse', 'hypre'
2244 
2245    Level: intermediate
2246 
2247 .seealso: `PCMGGalerkinSetMatProductAlgorithm()`
2248 
2249 @*/
2250 PetscErrorCode PCMGGalerkinGetMatProductAlgorithm(PC pc,const char *name[])
2251 {
2252   PetscFunctionBegin;
2253   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
2254   PetscTryMethod(pc,"PCMGGalerkinGetMatProductAlgorithm_C",(PC,const char*[]),(pc,name));
2255   PetscFunctionReturn(0);
2256 }
2257 
2258 /*MC
2259      PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre
2260 
2261    Options Database Keys:
2262 +   -pc_hypre_type - One of euclid, pilut, parasails, boomeramg, ams, ads
2263 .   -pc_hypre_boomeramg_nodal_coarsen <n> - where n is from 1 to 6 (see HYPRE_BOOMERAMGSetNodal())
2264 .   -pc_hypre_boomeramg_vec_interp_variant <v> - where v is from 1 to 3 (see HYPRE_BoomerAMGSetInterpVecVariant())
2265 -   Many others, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX preconditioner
2266 
2267    Level: intermediate
2268 
2269    Notes:
2270     Apart from pc_hypre_type (for which there is PCHYPRESetType()),
2271           the many hypre options can ONLY be set via the options database (e.g. the command line
2272           or with PetscOptionsSetValue(), there are no functions to set them)
2273 
2274           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_tol refer to the number of iterations
2275           (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if
2276           -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner
2277           (-pc_hypre_boomeramg_tol should be set to 0.0 - the default - to strictly use a fixed number of
2278           iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations
2279           and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10
2280           then AT MOST twenty V-cycles of boomeramg will be called.
2281 
2282            Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation
2283            (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry.
2284            Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi.
2285           If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly
2286           and use -ksp_max_it to control the number of V-cycles.
2287           (see the PETSc FAQ.html at the PETSc website under the Documentation tab).
2288 
2289           2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option
2290           -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L.
2291 
2292           MatSetNearNullSpace() - if you provide a near null space to your matrix it is ignored by hypre UNLESS you also use
2293           the following two options:
2294 
2295           See PCPFMG for access to the hypre Struct PFMG solver
2296 
2297    GPU Notes:
2298      To configure hypre BoomerAMG so that it can utilize NVIDIA GPUs run ./configure --download-hypre --with-cuda
2299      Then pass VECCUDA vectors and MATAIJCUSPARSE matrices to the solvers and PETSc will automatically utilize hypre's GPU solvers.
2300 
2301      To configure hypre BoomerAMG so that it can utilize AMD GPUs run ./configure --download-hypre --with-hip
2302      Then pass VECHIP vectors to the solvers and PETSc will automatically utilize hypre's GPU solvers.
2303 
2304 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
2305           `PCHYPRESetType()`, `PCPFMG`
2306 
2307 M*/
2308 
2309 PETSC_EXTERN PetscErrorCode PCCreate_HYPRE(PC pc)
2310 {
2311   PC_HYPRE       *jac;
2312 
2313   PetscFunctionBegin;
2314   PetscCall(PetscNewLog(pc,&jac));
2315 
2316   pc->data                = jac;
2317   pc->ops->reset          = PCReset_HYPRE;
2318   pc->ops->destroy        = PCDestroy_HYPRE;
2319   pc->ops->setfromoptions = PCSetFromOptions_HYPRE;
2320   pc->ops->setup          = PCSetUp_HYPRE;
2321   pc->ops->apply          = PCApply_HYPRE;
2322   jac->comm_hypre         = MPI_COMM_NULL;
2323   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetType_C",PCHYPRESetType_HYPRE));
2324   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPREGetType_C",PCHYPREGetType_HYPRE));
2325   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_HYPRE));
2326   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteGradient_C",PCHYPRESetDiscreteGradient_HYPRE));
2327   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteCurl_C",PCHYPRESetDiscreteCurl_HYPRE));
2328   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetInterpolations_C",PCHYPRESetInterpolations_HYPRE));
2329   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetEdgeConstantVectors_C",PCHYPRESetEdgeConstantVectors_HYPRE));
2330   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPREAMSSetInteriorNodes_C",PCHYPREAMSSetInteriorNodes_HYPRE));
2331   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetPoissonMatrix_C",PCHYPRESetPoissonMatrix_HYPRE));
2332   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGalerkinSetMatProductAlgorithm_C",PCMGGalerkinSetMatProductAlgorithm_HYPRE_BoomerAMG));
2333   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGalerkinGetMatProductAlgorithm_C",PCMGGalerkinGetMatProductAlgorithm_HYPRE_BoomerAMG));
2334 #if defined(PETSC_HAVE_HYPRE_DEVICE)
2335 #if defined(HYPRE_USING_HIP)
2336   PetscCall(PetscDeviceInitialize(PETSC_DEVICE_HIP));
2337 #endif
2338 #if defined(HYPRE_USING_CUDA)
2339   PetscCall(PetscDeviceInitialize(PETSC_DEVICE_CUDA));
2340 #endif
2341 #endif
2342   PetscFunctionReturn(0);
2343 }
2344 
2345 /* ---------------------------------------------------------------------------------------------------------------------------------*/
2346 
2347 typedef struct {
2348   MPI_Comm           hcomm;        /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */
2349   HYPRE_StructSolver hsolver;
2350 
2351   /* keep copy of PFMG options used so may view them */
2352   PetscInt its;
2353   double   tol;
2354   PetscInt relax_type;
2355   PetscInt rap_type;
2356   PetscInt num_pre_relax,num_post_relax;
2357   PetscInt max_levels;
2358   PetscInt skip_relax;
2359   PetscBool print_statistics;
2360 } PC_PFMG;
2361 
2362 PetscErrorCode PCDestroy_PFMG(PC pc)
2363 {
2364   PC_PFMG        *ex = (PC_PFMG*) pc->data;
2365 
2366   PetscFunctionBegin;
2367   if (ex->hsolver) PetscCallExternal(HYPRE_StructPFMGDestroy,ex->hsolver);
2368   PetscCall(PetscCommRestoreComm(PetscObjectComm((PetscObject)pc),&ex->hcomm));
2369   PetscCall(PetscFree(pc->data));
2370   PetscFunctionReturn(0);
2371 }
2372 
2373 static const char *PFMGRelaxType[] = {"Jacobi","Weighted-Jacobi","symmetric-Red/Black-Gauss-Seidel","Red/Black-Gauss-Seidel"};
2374 static const char *PFMGRAPType[] = {"Galerkin","non-Galerkin"};
2375 
2376 PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer)
2377 {
2378   PetscBool      iascii;
2379   PC_PFMG        *ex = (PC_PFMG*) pc->data;
2380 
2381   PetscFunctionBegin;
2382   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
2383   if (iascii) {
2384     PetscCall(PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG preconditioning\n"));
2385     PetscCall(PetscViewerASCIIPrintf(viewer,"    max iterations %" PetscInt_FMT "\n",ex->its));
2386     PetscCall(PetscViewerASCIIPrintf(viewer,"    tolerance %g\n",ex->tol));
2387     PetscCall(PetscViewerASCIIPrintf(viewer,"    relax type %s\n",PFMGRelaxType[ex->relax_type]));
2388     PetscCall(PetscViewerASCIIPrintf(viewer,"    RAP type %s\n",PFMGRAPType[ex->rap_type]));
2389     PetscCall(PetscViewerASCIIPrintf(viewer,"    number pre-relax %" PetscInt_FMT " post-relax %" PetscInt_FMT "\n",ex->num_pre_relax,ex->num_post_relax));
2390     PetscCall(PetscViewerASCIIPrintf(viewer,"    max levels %" PetscInt_FMT "\n",ex->max_levels));
2391     PetscCall(PetscViewerASCIIPrintf(viewer,"    skip relax %" PetscInt_FMT "\n",ex->skip_relax));
2392   }
2393   PetscFunctionReturn(0);
2394 }
2395 
2396 PetscErrorCode PCSetFromOptions_PFMG(PetscOptionItems *PetscOptionsObject,PC pc)
2397 {
2398   PC_PFMG        *ex = (PC_PFMG*) pc->data;
2399 
2400   PetscFunctionBegin;
2401   PetscOptionsHeadBegin(PetscOptionsObject,"PFMG options");
2402   PetscCall(PetscOptionsBool("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",ex->print_statistics,&ex->print_statistics,NULL));
2403   PetscCall(PetscOptionsInt("-pc_pfmg_its","Number of iterations of PFMG to use as preconditioner","HYPRE_StructPFMGSetMaxIter",ex->its,&ex->its,NULL));
2404   PetscCallExternal(HYPRE_StructPFMGSetMaxIter,ex->hsolver,ex->its);
2405   PetscCall(PetscOptionsInt("-pc_pfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_StructPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,NULL));
2406   PetscCallExternal(HYPRE_StructPFMGSetNumPreRelax,ex->hsolver,ex->num_pre_relax);
2407   PetscCall(PetscOptionsInt("-pc_pfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_StructPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,NULL));
2408   PetscCallExternal(HYPRE_StructPFMGSetNumPostRelax,ex->hsolver,ex->num_post_relax);
2409 
2410   PetscCall(PetscOptionsInt("-pc_pfmg_max_levels","Max Levels for MG hierarchy","HYPRE_StructPFMGSetMaxLevels",ex->max_levels,&ex->max_levels,NULL));
2411   PetscCallExternal(HYPRE_StructPFMGSetMaxLevels,ex->hsolver,ex->max_levels);
2412 
2413   PetscCall(PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,NULL));
2414   PetscCallExternal(HYPRE_StructPFMGSetTol,ex->hsolver,ex->tol);
2415   PetscCall(PetscOptionsEList("-pc_pfmg_relax_type","Relax type for the up and down cycles","HYPRE_StructPFMGSetRelaxType",PFMGRelaxType,PETSC_STATIC_ARRAY_LENGTH(PFMGRelaxType),PFMGRelaxType[ex->relax_type],&ex->relax_type,NULL));
2416   PetscCallExternal(HYPRE_StructPFMGSetRelaxType,ex->hsolver, ex->relax_type);
2417   PetscCall(PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,PETSC_STATIC_ARRAY_LENGTH(PFMGRAPType),PFMGRAPType[ex->rap_type],&ex->rap_type,NULL));
2418   PetscCallExternal(HYPRE_StructPFMGSetRAPType,ex->hsolver, ex->rap_type);
2419   PetscCall(PetscOptionsInt("-pc_pfmg_skip_relax","Skip relaxation on certain grids for isotropic problems. This can greatly improve efficiency by eliminating unnecessary relaxations when the underlying problem is isotropic","HYPRE_StructPFMGSetSkipRelax",ex->skip_relax,&ex->skip_relax,NULL));
2420   PetscCallExternal(HYPRE_StructPFMGSetSkipRelax,ex->hsolver, ex->skip_relax);
2421   PetscOptionsHeadEnd();
2422   PetscFunctionReturn(0);
2423 }
2424 
2425 PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y)
2426 {
2427   PC_PFMG           *ex = (PC_PFMG*) pc->data;
2428   PetscScalar       *yy;
2429   const PetscScalar *xx;
2430   PetscInt          ilower[3],iupper[3];
2431   HYPRE_Int         hlower[3],hupper[3];
2432   Mat_HYPREStruct   *mx = (Mat_HYPREStruct*)(pc->pmat->data);
2433 
2434   PetscFunctionBegin;
2435   PetscCall(PetscCitationsRegister(hypreCitation,&cite));
2436   PetscCall(DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]));
2437   /* when HYPRE_MIXEDINT is defined, sizeof(HYPRE_Int) == 32 */
2438   iupper[0] += ilower[0] - 1;
2439   iupper[1] += ilower[1] - 1;
2440   iupper[2] += ilower[2] - 1;
2441   hlower[0]  = (HYPRE_Int)ilower[0];
2442   hlower[1]  = (HYPRE_Int)ilower[1];
2443   hlower[2]  = (HYPRE_Int)ilower[2];
2444   hupper[0]  = (HYPRE_Int)iupper[0];
2445   hupper[1]  = (HYPRE_Int)iupper[1];
2446   hupper[2]  = (HYPRE_Int)iupper[2];
2447 
2448   /* copy x values over to hypre */
2449   PetscCallExternal(HYPRE_StructVectorSetConstantValues,mx->hb,0.0);
2450   PetscCall(VecGetArrayRead(x,&xx));
2451   PetscCallExternal(HYPRE_StructVectorSetBoxValues,mx->hb,hlower,hupper,(HYPRE_Complex*)xx);
2452   PetscCall(VecRestoreArrayRead(x,&xx));
2453   PetscCallExternal(HYPRE_StructVectorAssemble,mx->hb);
2454   PetscCallExternal(HYPRE_StructPFMGSolve,ex->hsolver,mx->hmat,mx->hb,mx->hx);
2455 
2456   /* copy solution values back to PETSc */
2457   PetscCall(VecGetArray(y,&yy));
2458   PetscCallExternal(HYPRE_StructVectorGetBoxValues,mx->hx,hlower,hupper,(HYPRE_Complex*)yy);
2459   PetscCall(VecRestoreArray(y,&yy));
2460   PetscFunctionReturn(0);
2461 }
2462 
2463 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)
2464 {
2465   PC_PFMG        *jac = (PC_PFMG*)pc->data;
2466   HYPRE_Int      oits;
2467 
2468   PetscFunctionBegin;
2469   PetscCall(PetscCitationsRegister(hypreCitation,&cite));
2470   PetscCallExternal(HYPRE_StructPFMGSetMaxIter,jac->hsolver,its*jac->its);
2471   PetscCallExternal(HYPRE_StructPFMGSetTol,jac->hsolver,rtol);
2472 
2473   PetscCall(PCApply_PFMG(pc,b,y));
2474   PetscCallExternal(HYPRE_StructPFMGGetNumIterations,jac->hsolver,&oits);
2475   *outits = oits;
2476   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
2477   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
2478   PetscCallExternal(HYPRE_StructPFMGSetTol,jac->hsolver,jac->tol);
2479   PetscCallExternal(HYPRE_StructPFMGSetMaxIter,jac->hsolver,jac->its);
2480   PetscFunctionReturn(0);
2481 }
2482 
2483 PetscErrorCode PCSetUp_PFMG(PC pc)
2484 {
2485   PC_PFMG         *ex = (PC_PFMG*) pc->data;
2486   Mat_HYPREStruct *mx = (Mat_HYPREStruct*)(pc->pmat->data);
2487   PetscBool       flg;
2488 
2489   PetscFunctionBegin;
2490   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESTRUCT,&flg));
2491   PetscCheck(flg,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESTRUCT with this preconditioner");
2492 
2493   /* create the hypre solver object and set its information */
2494   if (ex->hsolver) PetscCallExternal(HYPRE_StructPFMGDestroy,ex->hsolver);
2495   PetscCallExternal(HYPRE_StructPFMGCreate,ex->hcomm,&ex->hsolver);
2496 
2497   // Print Hypre statistics about the solve process
2498   if (ex->print_statistics) PetscCallExternal(HYPRE_StructPFMGSetPrintLevel,ex->hsolver,3);
2499 
2500   // The hypre options must be repeated here because the StructPFMG was destroyed and recreated
2501   PetscCallExternal(HYPRE_StructPFMGSetMaxIter     ,ex->hsolver, ex->its);
2502   PetscCallExternal(HYPRE_StructPFMGSetNumPreRelax ,ex->hsolver, ex->num_pre_relax);
2503   PetscCallExternal(HYPRE_StructPFMGSetNumPostRelax,ex->hsolver, ex->num_post_relax);
2504   PetscCallExternal(HYPRE_StructPFMGSetMaxLevels   ,ex->hsolver, ex->max_levels);
2505   PetscCallExternal(HYPRE_StructPFMGSetTol         ,ex->hsolver, ex->tol);
2506   PetscCallExternal(HYPRE_StructPFMGSetRelaxType   ,ex->hsolver, ex->relax_type);
2507   PetscCallExternal(HYPRE_StructPFMGSetRAPType     ,ex->hsolver, ex->rap_type);
2508 
2509   PetscCallExternal(HYPRE_StructPFMGSetup,ex->hsolver,mx->hmat,mx->hb,mx->hx);
2510   PetscCallExternal(HYPRE_StructPFMGSetZeroGuess,ex->hsolver);
2511   PetscFunctionReturn(0);
2512 }
2513 
2514 /*MC
2515      PCPFMG - the hypre PFMG multigrid solver
2516 
2517    Level: advanced
2518 
2519    Options Database:
2520 + -pc_pfmg_its <its> - number of iterations of PFMG to use as preconditioner
2521 . -pc_pfmg_num_pre_relax <steps> - number of smoothing steps before coarse grid solve
2522 . -pc_pfmg_num_post_relax <steps> - number of smoothing steps after coarse grid solve
2523 . -pc_pfmg_tol <tol> - tolerance of PFMG
2524 . -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
2525 . -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin
2526 - -pc_pfmg_skip_relax - skip relaxation on certain grids for isotropic problems. This can greatly improve efficiency by eliminating unnecessary relaxations when the underlying problem is isotropic, one of 0,1
2527 
2528    Notes:
2529     This is for CELL-centered descretizations
2530 
2531            This must be used with the MATHYPRESTRUCT matrix type.
2532            This is less general than in hypre, it supports only one block per process defined by a PETSc DMDA.
2533 
2534 .seealso: `PCMG`, `MATHYPRESTRUCT`
2535 M*/
2536 
2537 PETSC_EXTERN PetscErrorCode PCCreate_PFMG(PC pc)
2538 {
2539   PC_PFMG        *ex;
2540 
2541   PetscFunctionBegin;
2542   PetscCall(PetscNew(&ex)); \
2543   pc->data = ex;
2544 
2545   ex->its            = 1;
2546   ex->tol            = 1.e-8;
2547   ex->relax_type     = 1;
2548   ex->rap_type       = 0;
2549   ex->num_pre_relax  = 1;
2550   ex->num_post_relax = 1;
2551   ex->max_levels     = 0;
2552   ex->skip_relax     = 0;
2553   ex->print_statistics = PETSC_FALSE;
2554 
2555   pc->ops->setfromoptions  = PCSetFromOptions_PFMG;
2556   pc->ops->view            = PCView_PFMG;
2557   pc->ops->destroy         = PCDestroy_PFMG;
2558   pc->ops->apply           = PCApply_PFMG;
2559   pc->ops->applyrichardson = PCApplyRichardson_PFMG;
2560   pc->ops->setup           = PCSetUp_PFMG;
2561 
2562   PetscCall(PetscCommGetComm(PetscObjectComm((PetscObject)pc),&ex->hcomm));
2563   PetscCallExternal(HYPRE_StructPFMGCreate,ex->hcomm,&ex->hsolver);
2564   PetscFunctionReturn(0);
2565 }
2566 
2567 /* ---------------------------------------------------------------------------------------------------------------------------------------------------*/
2568 
2569 /* we know we are working with a HYPRE_SStructMatrix */
2570 typedef struct {
2571   MPI_Comm            hcomm;       /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */
2572   HYPRE_SStructSolver ss_solver;
2573 
2574   /* keep copy of SYSPFMG options used so may view them */
2575   PetscInt its;
2576   double   tol;
2577   PetscInt relax_type;
2578   PetscInt num_pre_relax,num_post_relax;
2579 } PC_SysPFMG;
2580 
2581 PetscErrorCode PCDestroy_SysPFMG(PC pc)
2582 {
2583   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;
2584 
2585   PetscFunctionBegin;
2586   if (ex->ss_solver) PetscCallExternal(HYPRE_SStructSysPFMGDestroy,ex->ss_solver);
2587   PetscCall(PetscCommRestoreComm(PetscObjectComm((PetscObject)pc),&ex->hcomm));
2588   PetscCall(PetscFree(pc->data));
2589   PetscFunctionReturn(0);
2590 }
2591 
2592 static const char *SysPFMGRelaxType[] = {"Weighted-Jacobi","Red/Black-Gauss-Seidel"};
2593 
2594 PetscErrorCode PCView_SysPFMG(PC pc,PetscViewer viewer)
2595 {
2596   PetscBool      iascii;
2597   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;
2598 
2599   PetscFunctionBegin;
2600   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
2601   if (iascii) {
2602     PetscCall(PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG preconditioning\n"));
2603     PetscCall(PetscViewerASCIIPrintf(viewer,"  max iterations %" PetscInt_FMT "\n",ex->its));
2604     PetscCall(PetscViewerASCIIPrintf(viewer,"  tolerance %g\n",ex->tol));
2605     PetscCall(PetscViewerASCIIPrintf(viewer,"  relax type %s\n",PFMGRelaxType[ex->relax_type]));
2606     PetscCall(PetscViewerASCIIPrintf(viewer,"  number pre-relax %" PetscInt_FMT " post-relax %" PetscInt_FMT "\n",ex->num_pre_relax,ex->num_post_relax));
2607   }
2608   PetscFunctionReturn(0);
2609 }
2610 
2611 PetscErrorCode PCSetFromOptions_SysPFMG(PetscOptionItems *PetscOptionsObject,PC pc)
2612 {
2613   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;
2614   PetscBool      flg = PETSC_FALSE;
2615 
2616   PetscFunctionBegin;
2617   PetscOptionsHeadBegin(PetscOptionsObject,"SysPFMG options");
2618   PetscCall(PetscOptionsBool("-pc_syspfmg_print_statistics","Print statistics","HYPRE_SStructSysPFMGSetPrintLevel",flg,&flg,NULL));
2619   if (flg) {
2620     PetscCallExternal(HYPRE_SStructSysPFMGSetPrintLevel,ex->ss_solver,3);
2621   }
2622   PetscCall(PetscOptionsInt("-pc_syspfmg_its","Number of iterations of SysPFMG to use as preconditioner","HYPRE_SStructSysPFMGSetMaxIter",ex->its,&ex->its,NULL));
2623   PetscCallExternal(HYPRE_SStructSysPFMGSetMaxIter,ex->ss_solver,ex->its);
2624   PetscCall(PetscOptionsInt("-pc_syspfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_SStructSysPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,NULL));
2625   PetscCallExternal(HYPRE_SStructSysPFMGSetNumPreRelax,ex->ss_solver,ex->num_pre_relax);
2626   PetscCall(PetscOptionsInt("-pc_syspfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_SStructSysPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,NULL));
2627   PetscCallExternal(HYPRE_SStructSysPFMGSetNumPostRelax,ex->ss_solver,ex->num_post_relax);
2628 
2629   PetscCall(PetscOptionsReal("-pc_syspfmg_tol","Tolerance of SysPFMG","HYPRE_SStructSysPFMGSetTol",ex->tol,&ex->tol,NULL));
2630   PetscCallExternal(HYPRE_SStructSysPFMGSetTol,ex->ss_solver,ex->tol);
2631   PetscCall(PetscOptionsEList("-pc_syspfmg_relax_type","Relax type for the up and down cycles","HYPRE_SStructSysPFMGSetRelaxType",SysPFMGRelaxType,PETSC_STATIC_ARRAY_LENGTH(SysPFMGRelaxType),SysPFMGRelaxType[ex->relax_type],&ex->relax_type,NULL));
2632   PetscCallExternal(HYPRE_SStructSysPFMGSetRelaxType,ex->ss_solver, ex->relax_type);
2633   PetscOptionsHeadEnd();
2634   PetscFunctionReturn(0);
2635 }
2636 
2637 PetscErrorCode PCApply_SysPFMG(PC pc,Vec x,Vec y)
2638 {
2639   PC_SysPFMG        *ex = (PC_SysPFMG*) pc->data;
2640   PetscScalar       *yy;
2641   const PetscScalar *xx;
2642   PetscInt          ilower[3],iupper[3];
2643   HYPRE_Int         hlower[3],hupper[3];
2644   Mat_HYPRESStruct  *mx     = (Mat_HYPRESStruct*)(pc->pmat->data);
2645   PetscInt          ordering= mx->dofs_order;
2646   PetscInt          nvars   = mx->nvars;
2647   PetscInt          part    = 0;
2648   PetscInt          size;
2649   PetscInt          i;
2650 
2651   PetscFunctionBegin;
2652   PetscCall(PetscCitationsRegister(hypreCitation,&cite));
2653   PetscCall(DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]));
2654   /* when HYPRE_MIXEDINT is defined, sizeof(HYPRE_Int) == 32 */
2655   iupper[0] += ilower[0] - 1;
2656   iupper[1] += ilower[1] - 1;
2657   iupper[2] += ilower[2] - 1;
2658   hlower[0]  = (HYPRE_Int)ilower[0];
2659   hlower[1]  = (HYPRE_Int)ilower[1];
2660   hlower[2]  = (HYPRE_Int)ilower[2];
2661   hupper[0]  = (HYPRE_Int)iupper[0];
2662   hupper[1]  = (HYPRE_Int)iupper[1];
2663   hupper[2]  = (HYPRE_Int)iupper[2];
2664 
2665   size = 1;
2666   for (i= 0; i< 3; i++) size *= (iupper[i]-ilower[i]+1);
2667 
2668   /* copy x values over to hypre for variable ordering */
2669   if (ordering) {
2670     PetscCallExternal(HYPRE_SStructVectorSetConstantValues,mx->ss_b,0.0);
2671     PetscCall(VecGetArrayRead(x,&xx));
2672     for (i= 0; i< nvars; i++) PetscCallExternal(HYPRE_SStructVectorSetBoxValues,mx->ss_b,part,hlower,hupper,i,(HYPRE_Complex*)(xx+(size*i)));
2673     PetscCall(VecRestoreArrayRead(x,&xx));
2674     PetscCallExternal(HYPRE_SStructVectorAssemble,mx->ss_b);
2675     PetscCallExternal(HYPRE_SStructMatrixMatvec,1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x);
2676     PetscCallExternal(HYPRE_SStructSysPFMGSolve,ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x);
2677 
2678     /* copy solution values back to PETSc */
2679     PetscCall(VecGetArray(y,&yy));
2680     for (i= 0; i< nvars; i++) PetscCallExternal(HYPRE_SStructVectorGetBoxValues,mx->ss_x,part,hlower,hupper,i,(HYPRE_Complex*)(yy+(size*i)));
2681     PetscCall(VecRestoreArray(y,&yy));
2682   } else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
2683     PetscScalar *z;
2684     PetscInt    j, k;
2685 
2686     PetscCall(PetscMalloc1(nvars*size,&z));
2687     PetscCallExternal(HYPRE_SStructVectorSetConstantValues,mx->ss_b,0.0);
2688     PetscCall(VecGetArrayRead(x,&xx));
2689 
2690     /* transform nodal to hypre's variable ordering for sys_pfmg */
2691     for (i= 0; i< size; i++) {
2692       k= i*nvars;
2693       for (j= 0; j< nvars; j++) z[j*size+i]= xx[k+j];
2694     }
2695     for (i= 0; i< nvars; i++) PetscCallExternal(HYPRE_SStructVectorSetBoxValues,mx->ss_b,part,hlower,hupper,i,(HYPRE_Complex*)(z+(size*i)));
2696     PetscCall(VecRestoreArrayRead(x,&xx));
2697     PetscCallExternal(HYPRE_SStructVectorAssemble,mx->ss_b);
2698     PetscCallExternal(HYPRE_SStructSysPFMGSolve,ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x);
2699 
2700     /* copy solution values back to PETSc */
2701     PetscCall(VecGetArray(y,&yy));
2702     for (i= 0; i< nvars; i++) PetscCallExternal(HYPRE_SStructVectorGetBoxValues,mx->ss_x,part,hlower,hupper,i,(HYPRE_Complex*)(z+(size*i)));
2703     /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
2704     for (i= 0; i< size; i++) {
2705       k= i*nvars;
2706       for (j= 0; j< nvars; j++) yy[k+j]= z[j*size+i];
2707     }
2708     PetscCall(VecRestoreArray(y,&yy));
2709     PetscCall(PetscFree(z));
2710   }
2711   PetscFunctionReturn(0);
2712 }
2713 
2714 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)
2715 {
2716   PC_SysPFMG     *jac = (PC_SysPFMG*)pc->data;
2717   HYPRE_Int      oits;
2718 
2719   PetscFunctionBegin;
2720   PetscCall(PetscCitationsRegister(hypreCitation,&cite));
2721   PetscCallExternal(HYPRE_SStructSysPFMGSetMaxIter,jac->ss_solver,its*jac->its);
2722   PetscCallExternal(HYPRE_SStructSysPFMGSetTol,jac->ss_solver,rtol);
2723   PetscCall(PCApply_SysPFMG(pc,b,y));
2724   PetscCallExternal(HYPRE_SStructSysPFMGGetNumIterations,jac->ss_solver,&oits);
2725   *outits = oits;
2726   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
2727   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
2728   PetscCallExternal(HYPRE_SStructSysPFMGSetTol,jac->ss_solver,jac->tol);
2729   PetscCallExternal(HYPRE_SStructSysPFMGSetMaxIter,jac->ss_solver,jac->its);
2730   PetscFunctionReturn(0);
2731 }
2732 
2733 PetscErrorCode PCSetUp_SysPFMG(PC pc)
2734 {
2735   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
2736   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct*)(pc->pmat->data);
2737   PetscBool        flg;
2738 
2739   PetscFunctionBegin;
2740   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESSTRUCT,&flg));
2741   PetscCheck(flg,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESSTRUCT with this preconditioner");
2742 
2743   /* create the hypre sstruct solver object and set its information */
2744   if (ex->ss_solver) PetscCallExternal(HYPRE_SStructSysPFMGDestroy,ex->ss_solver);
2745   PetscCallExternal(HYPRE_SStructSysPFMGCreate,ex->hcomm,&ex->ss_solver);
2746   PetscCallExternal(HYPRE_SStructSysPFMGSetZeroGuess,ex->ss_solver);
2747   PetscCallExternal(HYPRE_SStructSysPFMGSetup,ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x);
2748   PetscFunctionReturn(0);
2749 }
2750 
2751 /*MC
2752      PCSysPFMG - the hypre SysPFMG multigrid solver
2753 
2754    Level: advanced
2755 
2756    Options Database:
2757 + -pc_syspfmg_its <its> - number of iterations of SysPFMG to use as preconditioner
2758 . -pc_syspfmg_num_pre_relax <steps> - number of smoothing steps before coarse grid
2759 . -pc_syspfmg_num_post_relax <steps> - number of smoothing steps after coarse grid
2760 . -pc_syspfmg_tol <tol> - tolerance of SysPFMG
2761 - -pc_syspfmg_relax_type <Weighted-Jacobi,Red/Black-Gauss-Seidel> - relaxation type for the up and down cycles
2762 
2763    Notes:
2764     This is for CELL-centered descretizations
2765 
2766            This must be used with the MATHYPRESSTRUCT matrix type.
2767            This is less general than in hypre, it supports only one part, and one block per process defined by a PETSc DMDA.
2768            Also, only cell-centered variables.
2769 
2770 .seealso: `PCMG`, `MATHYPRESSTRUCT`
2771 M*/
2772 
2773 PETSC_EXTERN PetscErrorCode PCCreate_SysPFMG(PC pc)
2774 {
2775   PC_SysPFMG     *ex;
2776 
2777   PetscFunctionBegin;
2778   PetscCall(PetscNew(&ex)); \
2779   pc->data = ex;
2780 
2781   ex->its            = 1;
2782   ex->tol            = 1.e-8;
2783   ex->relax_type     = 1;
2784   ex->num_pre_relax  = 1;
2785   ex->num_post_relax = 1;
2786 
2787   pc->ops->setfromoptions  = PCSetFromOptions_SysPFMG;
2788   pc->ops->view            = PCView_SysPFMG;
2789   pc->ops->destroy         = PCDestroy_SysPFMG;
2790   pc->ops->apply           = PCApply_SysPFMG;
2791   pc->ops->applyrichardson = PCApplyRichardson_SysPFMG;
2792   pc->ops->setup           = PCSetUp_SysPFMG;
2793 
2794   PetscCall(PetscCommGetComm(PetscObjectComm((PetscObject)pc),&ex->hcomm));
2795   PetscCallExternal(HYPRE_SStructSysPFMGCreate,ex->hcomm,&ex->ss_solver);
2796   PetscFunctionReturn(0);
2797 }
2798 
2799 /* ---------------------------------------------------------------------------------------------------------------------------------*/
2800 
2801 // PC SMG
2802 typedef struct {
2803   MPI_Comm           hcomm;        /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */
2804   HYPRE_StructSolver hsolver;
2805   PetscInt           its;          /* keep copy of SMG options used so may view them */
2806   double             tol;
2807   PetscBool          print_statistics;
2808   PetscInt           num_pre_relax,num_post_relax;
2809 } PC_SMG;
2810 
2811 PetscErrorCode PCDestroy_SMG(PC pc)
2812 {
2813   PC_SMG        *ex = (PC_SMG*) pc->data;
2814 
2815   PetscFunctionBegin;
2816   if (ex->hsolver) PetscCallExternal(HYPRE_StructSMGDestroy,ex->hsolver);
2817   PetscCall(PetscCommRestoreComm(PetscObjectComm((PetscObject)pc),&ex->hcomm));
2818   PetscCall(PetscFree(pc->data));
2819   PetscFunctionReturn(0);
2820 }
2821 
2822 PetscErrorCode PCView_SMG(PC pc,PetscViewer viewer)
2823 {
2824   PetscBool      iascii;
2825   PC_SMG         *ex = (PC_SMG*) pc->data;
2826 
2827   PetscFunctionBegin;
2828   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
2829   if (iascii) {
2830     PetscCall(PetscViewerASCIIPrintf(viewer,"  HYPRE SMG preconditioning\n"));
2831     PetscCall(PetscViewerASCIIPrintf(viewer,"    max iterations %" PetscInt_FMT "\n",ex->its));
2832     PetscCall(PetscViewerASCIIPrintf(viewer,"    tolerance %g\n",ex->tol));
2833     PetscCall(PetscViewerASCIIPrintf(viewer,"    number pre-relax %" PetscInt_FMT " post-relax %" PetscInt_FMT "\n",ex->num_pre_relax,ex->num_post_relax));
2834   }
2835   PetscFunctionReturn(0);
2836 }
2837 
2838 PetscErrorCode PCSetFromOptions_SMG(PetscOptionItems *PetscOptionsObject,PC pc)
2839 {
2840   PC_SMG        *ex = (PC_SMG*) pc->data;
2841 
2842   PetscFunctionBegin;
2843   PetscOptionsHeadBegin(PetscOptionsObject,"SMG options");
2844 
2845   PetscCall(PetscOptionsInt("-pc_smg_its","Number of iterations of SMG to use as preconditioner","HYPRE_StructSMGSetMaxIter",ex->its,&ex->its,NULL));
2846   PetscCall(PetscOptionsInt("-pc_smg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_StructSMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,NULL));
2847   PetscCall(PetscOptionsInt("-pc_smg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_StructSMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,NULL));
2848   PetscCall(PetscOptionsReal("-pc_smg_tol","Tolerance of SMG","HYPRE_StructSMGSetTol",ex->tol,&ex->tol,NULL));
2849 
2850   PetscOptionsHeadEnd();
2851   PetscFunctionReturn(0);
2852 }
2853 
2854 PetscErrorCode PCApply_SMG(PC pc,Vec x,Vec y)
2855 {
2856   PC_SMG           *ex = (PC_SMG*) pc->data;
2857   PetscScalar       *yy;
2858   const PetscScalar *xx;
2859   PetscInt          ilower[3],iupper[3];
2860   HYPRE_Int         hlower[3],hupper[3];
2861   Mat_HYPREStruct   *mx = (Mat_HYPREStruct*)(pc->pmat->data);
2862 
2863   PetscFunctionBegin;
2864   PetscCall(PetscCitationsRegister(hypreCitation,&cite));
2865   PetscCall(DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]));
2866   /* when HYPRE_MIXEDINT is defined, sizeof(HYPRE_Int) == 32 */
2867   iupper[0] += ilower[0] - 1;
2868   iupper[1] += ilower[1] - 1;
2869   iupper[2] += ilower[2] - 1;
2870   hlower[0]  = (HYPRE_Int)ilower[0];
2871   hlower[1]  = (HYPRE_Int)ilower[1];
2872   hlower[2]  = (HYPRE_Int)ilower[2];
2873   hupper[0]  = (HYPRE_Int)iupper[0];
2874   hupper[1]  = (HYPRE_Int)iupper[1];
2875   hupper[2]  = (HYPRE_Int)iupper[2];
2876 
2877   /* copy x values over to hypre */
2878   PetscCallExternal(HYPRE_StructVectorSetConstantValues,mx->hb,0.0);
2879   PetscCall(VecGetArrayRead(x,&xx));
2880   PetscCallExternal(HYPRE_StructVectorSetBoxValues,mx->hb,hlower,hupper,(HYPRE_Complex*)xx);
2881   PetscCall(VecRestoreArrayRead(x,&xx));
2882   PetscCallExternal(HYPRE_StructVectorAssemble,mx->hb);
2883   PetscCallExternal(HYPRE_StructSMGSolve,ex->hsolver,mx->hmat,mx->hb,mx->hx);
2884 
2885   /* copy solution values back to PETSc */
2886   PetscCall(VecGetArray(y,&yy));
2887   PetscCallExternal(HYPRE_StructVectorGetBoxValues,mx->hx,hlower,hupper,(HYPRE_Complex*)yy);
2888   PetscCall(VecRestoreArray(y,&yy));
2889   PetscFunctionReturn(0);
2890 }
2891 
2892 static PetscErrorCode PCApplyRichardson_SMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
2893 {
2894   PC_SMG        *jac = (PC_SMG*)pc->data;
2895   HYPRE_Int      oits;
2896 
2897   PetscFunctionBegin;
2898   PetscCall(PetscCitationsRegister(hypreCitation,&cite));
2899   PetscCallExternal(HYPRE_StructSMGSetMaxIter,jac->hsolver,its*jac->its);
2900   PetscCallExternal(HYPRE_StructSMGSetTol,jac->hsolver,rtol);
2901 
2902   PetscCall(PCApply_SMG(pc,b,y));
2903   PetscCallExternal(HYPRE_StructSMGGetNumIterations,jac->hsolver,&oits);
2904   *outits = oits;
2905   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
2906   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
2907   PetscCallExternal(HYPRE_StructSMGSetTol,jac->hsolver,jac->tol);
2908   PetscCallExternal(HYPRE_StructSMGSetMaxIter,jac->hsolver,jac->its);
2909   PetscFunctionReturn(0);
2910 }
2911 
2912 PetscErrorCode PCSetUp_SMG(PC pc)
2913 {
2914   PetscInt        i, dim;
2915   PC_SMG         *ex = (PC_SMG*) pc->data;
2916   Mat_HYPREStruct *mx = (Mat_HYPREStruct*)(pc->pmat->data);
2917   PetscBool       flg;
2918   DMBoundaryType  p[3];
2919   PetscInt        M[3];
2920 
2921   PetscFunctionBegin;
2922   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESTRUCT,&flg));
2923   PetscCheck(flg,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESTRUCT with this preconditioner");
2924 
2925   PetscCall(DMDAGetInfo(mx->da,&dim,&M[0], &M[1], &M[2],0,0,0,0,0,&p[0],&p[1],&p[2],0));
2926   // Check if power of 2 in periodic directions
2927   for (i=0;i<dim;i++){
2928     if (((M[i] & (M[i] - 1)) != 0) && (p[i]==DM_BOUNDARY_PERIODIC)) {
2929       SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"With SMG, the number of points in a periodic direction must be a power of 2, but is here %" PetscInt_FMT ".",M[i]);
2930     }
2931   }
2932 
2933   /* create the hypre solver object and set its information */
2934   if (ex->hsolver) PetscCallExternal(HYPRE_StructSMGDestroy,(ex->hsolver));
2935   PetscCallExternal(HYPRE_StructSMGCreate,ex->hcomm,&ex->hsolver);
2936   // The hypre options must be set here and not in SetFromOptions because it is created here!
2937   PetscCallExternal(HYPRE_StructSMGSetMaxIter,ex->hsolver,ex->its);
2938   PetscCallExternal(HYPRE_StructSMGSetNumPreRelax,ex->hsolver,ex->num_pre_relax);
2939   PetscCallExternal(HYPRE_StructSMGSetNumPostRelax,ex->hsolver,ex->num_post_relax);
2940   PetscCallExternal(HYPRE_StructSMGSetTol,ex->hsolver,ex->tol);
2941 
2942   PetscCallExternal(HYPRE_StructSMGSetup,ex->hsolver,mx->hmat,mx->hb,mx->hx);
2943   PetscCallExternal(HYPRE_StructSMGSetZeroGuess,ex->hsolver);
2944   PetscFunctionReturn(0);
2945 }
2946 
2947 /*MC
2948      PCSMG - the hypre SMG multigrid solver
2949 
2950    Level: advanced
2951 
2952    Options Database:
2953 + -pc_smg_its <its> number of iterations of SMG to use as preconditioner
2954 . -pc_smg_num_pre_relax <steps> number of smoothing steps before coarse grid
2955 . -pc_smg_num_post_relax <steps> number of smoothing steps after coarse grid
2956 . -pc_smg_tol <tol> tolerance of SMG
2957 
2958    Notes:
2959     This is for CELL-centered descretizations
2960 
2961            This must be used with the MATHYPRESTRUCT matrix type.
2962            This is less general than in hypre, it supports only one block per process defined by a PETSc DMDA.
2963 
2964 .seealso:  PCMG, MATHYPRESTRUCT, PCPFMG
2965 M*/
2966 
2967 PETSC_EXTERN PetscErrorCode PCCreate_SMG(PC pc)
2968 {
2969   PC_SMG        *ex;
2970 
2971   PetscFunctionBegin;
2972   PetscCall(PetscNew(&ex)); \
2973   pc->data = ex;
2974 
2975   ex->its            = 1;
2976   ex->tol            = 1.e-8;
2977   ex->num_pre_relax  = 1;
2978   ex->num_post_relax = 1;
2979 
2980   pc->ops->setfromoptions  = PCSetFromOptions_SMG;
2981   pc->ops->view            = PCView_SMG;
2982   pc->ops->destroy         = PCDestroy_SMG;
2983   pc->ops->apply           = PCApply_SMG;
2984   pc->ops->applyrichardson = PCApplyRichardson_SMG;
2985   pc->ops->setup           = PCSetUp_SMG;
2986 
2987   PetscCall(PetscCommGetComm(PetscObjectComm((PetscObject)pc),&ex->hcomm));
2988   PetscCallExternal(HYPRE_StructSMGCreate,ex->hcomm,&ex->hsolver);
2989   PetscFunctionReturn(0);
2990 }
2991