xref: /petsc/src/ksp/pc/impls/hypre/hypre.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
116d9e3a6SLisandro Dalcin /*
216d9e3a6SLisandro Dalcin    Provides an interface to the LLNL package hypre
316d9e3a6SLisandro Dalcin */
40f1074feSSatish Balay 
5589dcaf0SStefano Zampini #include <petscpkg_version.h>
6af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
749a781f5SStefano Zampini /* this include is needed ONLY to allow access to the private data inside the Mat object specific to hypre */
849a781f5SStefano Zampini #include <petsc/private/matimpl.h>
96ea7df73SStefano Zampini #include <petsc/private/vecimpl.h>
1058968eb6SStefano Zampini #include <../src/vec/vec/impls/hypre/vhyp.h>
1149a781f5SStefano Zampini #include <../src/mat/impls/hypre/mhypre.h>
12c6db04a5SJed Brown #include <../src/dm/impls/da/hypre/mhyp.h>
134cb006feSStefano Zampini #include <_hypre_parcsr_ls.h>
148a2c336bSFande Kong #include <petscmathypre.h>
1516d9e3a6SLisandro Dalcin 
16a4af0ceeSJacob Faibussowitsch #if defined(PETSC_HAVE_HYPRE_DEVICE)
17a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h>
18a4af0ceeSJacob Faibussowitsch #endif
19a4af0ceeSJacob Faibussowitsch 
20dff31646SBarry Smith static PetscBool  cite            = PETSC_FALSE;
219371c9d4SSatish Balay static const char hypreCitation[] = "@manual{hypre-web-page,\n  title  = {{\\sl hypre}: High Performance Preconditioners},\n  organization = {Lawrence Livermore National Laboratory},\n  note  = "
229371c9d4SSatish Balay                                     "{\\url{https://computation.llnl.gov/projects/hypre-scalable-linear-solvers-multigrid-methods}}\n}\n";
231f817a21SBarry Smith 
2416d9e3a6SLisandro Dalcin /*
2516d9e3a6SLisandro Dalcin    Private context (data structure) for the  preconditioner.
2616d9e3a6SLisandro Dalcin */
2716d9e3a6SLisandro Dalcin typedef struct {
2816d9e3a6SLisandro Dalcin   HYPRE_Solver hsolver;
2949a781f5SStefano Zampini   Mat          hpmat; /* MatHYPRE */
3016d9e3a6SLisandro Dalcin 
314ddd07fcSJed Brown   HYPRE_Int (*destroy)(HYPRE_Solver);
324ddd07fcSJed Brown   HYPRE_Int (*solve)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
334ddd07fcSJed Brown   HYPRE_Int (*setup)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
3416d9e3a6SLisandro Dalcin 
3516d9e3a6SLisandro Dalcin   MPI_Comm comm_hypre;
3616d9e3a6SLisandro Dalcin   char    *hypre_type;
3716d9e3a6SLisandro Dalcin 
3816d9e3a6SLisandro Dalcin   /* options for Pilut and BoomerAMG*/
394ddd07fcSJed Brown   PetscInt  maxiter;
4039accc25SStefano Zampini   PetscReal tol;
4116d9e3a6SLisandro Dalcin 
4216d9e3a6SLisandro Dalcin   /* options for Pilut */
434ddd07fcSJed Brown   PetscInt factorrowsize;
4416d9e3a6SLisandro Dalcin 
4516d9e3a6SLisandro Dalcin   /* options for ParaSails */
464ddd07fcSJed Brown   PetscInt  nlevels;
478966356dSPierre Jolivet   PetscReal threshold;
4839accc25SStefano Zampini   PetscReal filter;
4939accc25SStefano Zampini   PetscReal loadbal;
504ddd07fcSJed Brown   PetscInt  logging;
514ddd07fcSJed Brown   PetscInt  ruse;
524ddd07fcSJed Brown   PetscInt  symt;
5316d9e3a6SLisandro Dalcin 
5422b6d1caSBarry Smith   /* options for BoomerAMG */
55ace3abfcSBarry Smith   PetscBool printstatistics;
5616d9e3a6SLisandro Dalcin 
5716d9e3a6SLisandro Dalcin   /* options for BoomerAMG */
584ddd07fcSJed Brown   PetscInt  cycletype;
594ddd07fcSJed Brown   PetscInt  maxlevels;
6039accc25SStefano Zampini   PetscReal strongthreshold;
6139accc25SStefano Zampini   PetscReal maxrowsum;
624ddd07fcSJed Brown   PetscInt  gridsweeps[3];
634ddd07fcSJed Brown   PetscInt  coarsentype;
644ddd07fcSJed Brown   PetscInt  measuretype;
656a251517SEike Mueller   PetscInt  smoothtype;
668131ecf7SEike Mueller   PetscInt  smoothnumlevels;
67ec64516dSEike Mueller   PetscInt  eu_level;         /* Number of levels for ILU(k) in Euclid */
6839accc25SStefano Zampini   PetscReal eu_droptolerance; /* Drop tolerance for ILU(k) in Euclid */
69ec64516dSEike Mueller   PetscInt  eu_bj;            /* Defines use of Block Jacobi ILU in Euclid */
704ddd07fcSJed Brown   PetscInt  relaxtype[3];
7139accc25SStefano Zampini   PetscReal relaxweight;
7239accc25SStefano Zampini   PetscReal outerrelaxweight;
734ddd07fcSJed Brown   PetscInt  relaxorder;
7439accc25SStefano Zampini   PetscReal truncfactor;
75ace3abfcSBarry Smith   PetscBool applyrichardson;
764ddd07fcSJed Brown   PetscInt  pmax;
774ddd07fcSJed Brown   PetscInt  interptype;
78589dcaf0SStefano Zampini   PetscInt  maxc;
79589dcaf0SStefano Zampini   PetscInt  minc;
80db6f9c32SMark Adams #if PETSC_PKG_HYPRE_VERSION_GE(2, 23, 0)
81db6f9c32SMark Adams   char *spgemm_type; // this is a global hypre parameter but is closely associated with BoomerAMG
82db6f9c32SMark Adams #endif
836ea7df73SStefano Zampini   /* GPU */
846ea7df73SStefano Zampini   PetscBool keeptranspose;
856ea7df73SStefano Zampini   PetscInt  rap2;
866ea7df73SStefano Zampini   PetscInt  mod_rap2;
876ea7df73SStefano Zampini 
88589dcaf0SStefano Zampini   /* AIR */
89589dcaf0SStefano Zampini   PetscInt  Rtype;
90589dcaf0SStefano Zampini   PetscReal Rstrongthreshold;
91589dcaf0SStefano Zampini   PetscReal Rfilterthreshold;
92589dcaf0SStefano Zampini   PetscInt  Adroptype;
93589dcaf0SStefano Zampini   PetscReal Adroptol;
94589dcaf0SStefano Zampini 
954ddd07fcSJed Brown   PetscInt  agg_nl;
966ea7df73SStefano Zampini   PetscInt  agg_interptype;
974ddd07fcSJed Brown   PetscInt  agg_num_paths;
98ace3abfcSBarry Smith   PetscBool nodal_relax;
994ddd07fcSJed Brown   PetscInt  nodal_relax_levels;
1004cb006feSStefano Zampini 
1015272c319SBarry Smith   PetscInt  nodal_coarsening;
10222e51d31SStefano Zampini   PetscInt  nodal_coarsening_diag;
1035272c319SBarry Smith   PetscInt  vec_interp_variant;
10422e51d31SStefano Zampini   PetscInt  vec_interp_qmax;
10522e51d31SStefano Zampini   PetscBool vec_interp_smooth;
10622e51d31SStefano Zampini   PetscInt  interp_refine;
10722e51d31SStefano Zampini 
1086ea7df73SStefano Zampini   /* NearNullSpace support */
1096ea7df73SStefano Zampini   VecHYPRE_IJVector *hmnull;
1106ea7df73SStefano Zampini   HYPRE_ParVector   *phmnull;
1115272c319SBarry Smith   PetscInt           n_hmnull;
1125272c319SBarry Smith   Vec                hmnull_constant;
1135272c319SBarry Smith 
114863406b8SStefano Zampini   /* options for AS (Auxiliary Space preconditioners) */
115863406b8SStefano Zampini   PetscInt  as_print;
116863406b8SStefano Zampini   PetscInt  as_max_iter;
117863406b8SStefano Zampini   PetscReal as_tol;
118863406b8SStefano Zampini   PetscInt  as_relax_type;
119863406b8SStefano Zampini   PetscInt  as_relax_times;
120863406b8SStefano Zampini   PetscReal as_relax_weight;
121863406b8SStefano Zampini   PetscReal as_omega;
122863406b8SStefano Zampini   PetscInt  as_amg_alpha_opts[5]; /* AMG coarsen type, agg_levels, relax_type, interp_type, Pmax for vector Poisson (AMS) or Curl problem (ADS) */
123863406b8SStefano Zampini   PetscReal as_amg_alpha_theta;   /* AMG strength for vector Poisson (AMS) or Curl problem (ADS) */
124863406b8SStefano Zampini   PetscInt  as_amg_beta_opts[5];  /* AMG coarsen type, agg_levels, relax_type, interp_type, Pmax for scalar Poisson (AMS) or vector Poisson (ADS) */
125863406b8SStefano Zampini   PetscReal as_amg_beta_theta;    /* AMG strength for scalar Poisson (AMS) or vector Poisson (ADS)  */
1264cb006feSStefano Zampini   PetscInt  ams_cycle_type;
127863406b8SStefano Zampini   PetscInt  ads_cycle_type;
1284cb006feSStefano Zampini 
1294cb006feSStefano Zampini   /* additional data */
1305ac14e1cSStefano Zampini   Mat G;             /* MatHYPRE */
1315ac14e1cSStefano Zampini   Mat C;             /* MatHYPRE */
1325ac14e1cSStefano Zampini   Mat alpha_Poisson; /* MatHYPRE */
1335ac14e1cSStefano Zampini   Mat beta_Poisson;  /* MatHYPRE */
1345ac14e1cSStefano Zampini 
1355ac14e1cSStefano Zampini   /* extra information for AMS */
1365ac14e1cSStefano Zampini   PetscInt          dim; /* geometrical dimension */
1376ea7df73SStefano Zampini   VecHYPRE_IJVector coords[3];
1386ea7df73SStefano Zampini   VecHYPRE_IJVector constants[3];
139be14dc20SKerry Key   VecHYPRE_IJVector interior;
1406bf688a0SCe Qin   Mat               RT_PiFull, RT_Pi[3];
1416bf688a0SCe Qin   Mat               ND_PiFull, ND_Pi[3];
1424cb006feSStefano Zampini   PetscBool         ams_beta_is_zero;
14323df4f25SStefano Zampini   PetscBool         ams_beta_is_zero_part;
14423df4f25SStefano Zampini   PetscInt          ams_proj_freq;
14516d9e3a6SLisandro Dalcin } PC_HYPRE;
14616d9e3a6SLisandro Dalcin 
1479371c9d4SSatish Balay PetscErrorCode PCHYPREGetSolver(PC pc, HYPRE_Solver *hsolver) {
148d2128fa2SBarry Smith   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
149d2128fa2SBarry Smith 
150d2128fa2SBarry Smith   PetscFunctionBegin;
151d2128fa2SBarry Smith   *hsolver = jac->hsolver;
152d2128fa2SBarry Smith   PetscFunctionReturn(0);
153d2128fa2SBarry Smith }
15416d9e3a6SLisandro Dalcin 
155fd2dd295SFande Kong /*
1568a2c336bSFande Kong   Matrices with AIJ format are created IN PLACE with using (I,J,data) from BoomerAMG. Since the data format in hypre_ParCSRMatrix
1578a2c336bSFande Kong   is different from that used in PETSc, the original hypre_ParCSRMatrix can not be used any more after call this routine.
1588a2c336bSFande Kong   It is used in PCHMG. Other users should avoid using this function.
159fd2dd295SFande Kong */
1609371c9d4SSatish Balay static PetscErrorCode PCGetCoarseOperators_BoomerAMG(PC pc, PetscInt *nlevels, Mat *operators[]) {
1618a2c336bSFande Kong   PC_HYPRE            *jac  = (PC_HYPRE *)pc->data;
1628a2c336bSFande Kong   PetscBool            same = PETSC_FALSE;
1638a2c336bSFande Kong   PetscInt             num_levels, l;
1648a2c336bSFande Kong   Mat                 *mattmp;
1658a2c336bSFande Kong   hypre_ParCSRMatrix **A_array;
1668a2c336bSFande Kong 
1678a2c336bSFande Kong   PetscFunctionBegin;
1689566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(jac->hypre_type, "boomeramg", &same));
1695f80ce2aSJacob Faibussowitsch   PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_NOTSAMETYPE, "Hypre type is not BoomerAMG ");
1708a2c336bSFande Kong   num_levels = hypre_ParAMGDataNumLevels((hypre_ParAMGData *)(jac->hsolver));
1719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(num_levels, &mattmp));
1728a2c336bSFande Kong   A_array = hypre_ParAMGDataAArray((hypre_ParAMGData *)(jac->hsolver));
1738a2c336bSFande Kong   for (l = 1; l < num_levels; l++) {
1749566063dSJacob Faibussowitsch     PetscCall(MatCreateFromParCSR(A_array[l], MATAIJ, PETSC_OWN_POINTER, &(mattmp[num_levels - 1 - l])));
1758a2c336bSFande Kong     /* We want to own the data, and HYPRE can not touch this matrix any more */
1768a2c336bSFande Kong     A_array[l] = NULL;
1778a2c336bSFande Kong   }
1788a2c336bSFande Kong   *nlevels   = num_levels;
1798a2c336bSFande Kong   *operators = mattmp;
1808a2c336bSFande Kong   PetscFunctionReturn(0);
1818a2c336bSFande Kong }
1828a2c336bSFande Kong 
183fd2dd295SFande Kong /*
1848a2c336bSFande Kong   Matrices with AIJ format are created IN PLACE with using (I,J,data) from BoomerAMG. Since the data format in hypre_ParCSRMatrix
1858a2c336bSFande Kong   is different from that used in PETSc, the original hypre_ParCSRMatrix can not be used any more after call this routine.
1868a2c336bSFande Kong   It is used in PCHMG. Other users should avoid using this function.
187fd2dd295SFande Kong */
1889371c9d4SSatish Balay static PetscErrorCode PCGetInterpolations_BoomerAMG(PC pc, PetscInt *nlevels, Mat *interpolations[]) {
1898a2c336bSFande Kong   PC_HYPRE            *jac  = (PC_HYPRE *)pc->data;
1908a2c336bSFande Kong   PetscBool            same = PETSC_FALSE;
1918a2c336bSFande Kong   PetscInt             num_levels, l;
1928a2c336bSFande Kong   Mat                 *mattmp;
1938a2c336bSFande Kong   hypre_ParCSRMatrix **P_array;
1948a2c336bSFande Kong 
1958a2c336bSFande Kong   PetscFunctionBegin;
1969566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(jac->hypre_type, "boomeramg", &same));
1975f80ce2aSJacob Faibussowitsch   PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_NOTSAMETYPE, "Hypre type is not BoomerAMG ");
1988a2c336bSFande Kong   num_levels = hypre_ParAMGDataNumLevels((hypre_ParAMGData *)(jac->hsolver));
1999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(num_levels, &mattmp));
2008a2c336bSFande Kong   P_array = hypre_ParAMGDataPArray((hypre_ParAMGData *)(jac->hsolver));
2018a2c336bSFande Kong   for (l = 1; l < num_levels; l++) {
2029566063dSJacob Faibussowitsch     PetscCall(MatCreateFromParCSR(P_array[num_levels - 1 - l], MATAIJ, PETSC_OWN_POINTER, &(mattmp[l - 1])));
2038a2c336bSFande Kong     /* We want to own the data, and HYPRE can not touch this matrix any more */
2048a2c336bSFande Kong     P_array[num_levels - 1 - l] = NULL;
2058a2c336bSFande Kong   }
2068a2c336bSFande Kong   *nlevels        = num_levels;
2078a2c336bSFande Kong   *interpolations = mattmp;
2088a2c336bSFande Kong   PetscFunctionReturn(0);
2098a2c336bSFande Kong }
2108a2c336bSFande Kong 
211ce6a8a0dSJed Brown /* Resets (frees) Hypre's representation of the near null space */
2129371c9d4SSatish Balay static PetscErrorCode PCHYPREResetNearNullSpace_Private(PC pc) {
213ce6a8a0dSJed Brown   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
214ce6a8a0dSJed Brown   PetscInt  i;
215ce6a8a0dSJed Brown 
2169d678128SJed Brown   PetscFunctionBegin;
21748a46eb9SPierre Jolivet   for (i = 0; i < jac->n_hmnull; i++) PetscCall(VecHYPRE_IJVectorDestroy(&jac->hmnull[i]));
2189566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->hmnull));
2199566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->phmnull));
2209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->hmnull_constant));
2219d678128SJed Brown   jac->n_hmnull = 0;
222ce6a8a0dSJed Brown   PetscFunctionReturn(0);
223ce6a8a0dSJed Brown }
224ce6a8a0dSJed Brown 
2259371c9d4SSatish Balay static PetscErrorCode PCSetUp_HYPRE(PC pc) {
22616d9e3a6SLisandro Dalcin   PC_HYPRE          *jac = (PC_HYPRE *)pc->data;
22749a781f5SStefano Zampini   Mat_HYPRE         *hjac;
22816d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
22916d9e3a6SLisandro Dalcin   HYPRE_ParVector    bv, xv;
23049a781f5SStefano Zampini   PetscBool          ishypre;
23116d9e3a6SLisandro Dalcin 
23216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
23348a46eb9SPierre Jolivet   if (!jac->hypre_type) PetscCall(PCHYPRESetType(pc, "boomeramg"));
2345f5c5b43SBarry Smith 
2359566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATHYPRE, &ishypre));
23649a781f5SStefano Zampini   if (!ishypre) {
2379566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&jac->hpmat));
2389566063dSJacob Faibussowitsch     PetscCall(MatConvert(pc->pmat, MATHYPRE, MAT_INITIAL_MATRIX, &jac->hpmat));
23949a781f5SStefano Zampini   } else {
2409566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)pc->pmat));
2419566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&jac->hpmat));
24249a781f5SStefano Zampini     jac->hpmat = pc->pmat;
24316d9e3a6SLisandro Dalcin   }
2446ea7df73SStefano Zampini   /* allow debug */
2459566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(jac->hpmat, NULL, "-pc_hypre_mat_view"));
24649a781f5SStefano Zampini   hjac = (Mat_HYPRE *)(jac->hpmat->data);
2475f5c5b43SBarry Smith 
24816d9e3a6SLisandro Dalcin   /* special case for BoomerAMG */
24916d9e3a6SLisandro Dalcin   if (jac->setup == HYPRE_BoomerAMGSetup) {
2505272c319SBarry Smith     MatNullSpace mnull;
2515272c319SBarry Smith     PetscBool    has_const;
25249a781f5SStefano Zampini     PetscInt     bs, nvec, i;
2535272c319SBarry Smith     const Vec   *vecs;
2545272c319SBarry Smith 
2559566063dSJacob Faibussowitsch     PetscCall(MatGetBlockSize(pc->pmat, &bs));
256792fecdfSBarry Smith     if (bs > 1) PetscCallExternal(HYPRE_BoomerAMGSetNumFunctions, jac->hsolver, bs);
2579566063dSJacob Faibussowitsch     PetscCall(MatGetNearNullSpace(pc->mat, &mnull));
2585272c319SBarry Smith     if (mnull) {
2599566063dSJacob Faibussowitsch       PetscCall(PCHYPREResetNearNullSpace_Private(pc));
2609566063dSJacob Faibussowitsch       PetscCall(MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs));
2619566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nvec + 1, &jac->hmnull));
2629566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nvec + 1, &jac->phmnull));
2635272c319SBarry Smith       for (i = 0; i < nvec; i++) {
2649566063dSJacob Faibussowitsch         PetscCall(VecHYPRE_IJVectorCreate(vecs[i]->map, &jac->hmnull[i]));
2659566063dSJacob Faibussowitsch         PetscCall(VecHYPRE_IJVectorCopy(vecs[i], jac->hmnull[i]));
266792fecdfSBarry Smith         PetscCallExternal(HYPRE_IJVectorGetObject, jac->hmnull[i]->ij, (void **)&jac->phmnull[i]);
2675272c319SBarry Smith       }
2685272c319SBarry Smith       if (has_const) {
2699566063dSJacob Faibussowitsch         PetscCall(MatCreateVecs(pc->pmat, &jac->hmnull_constant, NULL));
2709566063dSJacob Faibussowitsch         PetscCall(VecSet(jac->hmnull_constant, 1));
2719566063dSJacob Faibussowitsch         PetscCall(VecNormalize(jac->hmnull_constant, NULL));
2729566063dSJacob Faibussowitsch         PetscCall(VecHYPRE_IJVectorCreate(jac->hmnull_constant->map, &jac->hmnull[nvec]));
2739566063dSJacob Faibussowitsch         PetscCall(VecHYPRE_IJVectorCopy(jac->hmnull_constant, jac->hmnull[nvec]));
274792fecdfSBarry Smith         PetscCallExternal(HYPRE_IJVectorGetObject, jac->hmnull[nvec]->ij, (void **)&jac->phmnull[nvec]);
2755272c319SBarry Smith         nvec++;
2765272c319SBarry Smith       }
277792fecdfSBarry Smith       PetscCallExternal(HYPRE_BoomerAMGSetInterpVectors, jac->hsolver, nvec, jac->phmnull);
2785272c319SBarry Smith       jac->n_hmnull = nvec;
2795272c319SBarry Smith     }
2804cb006feSStefano Zampini   }
281863406b8SStefano Zampini 
2824cb006feSStefano Zampini   /* special case for AMS */
2834cb006feSStefano Zampini   if (jac->setup == HYPRE_AMSSetup) {
2845ac14e1cSStefano Zampini     Mat_HYPRE         *hm;
2855ac14e1cSStefano Zampini     HYPRE_ParCSRMatrix parcsr;
2866bf688a0SCe Qin     if (!jac->coords[0] && !jac->constants[0] && !(jac->ND_PiFull || (jac->ND_Pi[0] && jac->ND_Pi[1]))) {
287f1580f4eSBarry Smith       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()");
2886bf688a0SCe Qin     }
28948a46eb9SPierre Jolivet     if (jac->dim) PetscCallExternal(HYPRE_AMSSetDimension, jac->hsolver, jac->dim);
2905ac14e1cSStefano Zampini     if (jac->constants[0]) {
2915ac14e1cSStefano Zampini       HYPRE_ParVector ozz, zoz, zzo = NULL;
292792fecdfSBarry Smith       PetscCallExternal(HYPRE_IJVectorGetObject, jac->constants[0]->ij, (void **)(&ozz));
293792fecdfSBarry Smith       PetscCallExternal(HYPRE_IJVectorGetObject, jac->constants[1]->ij, (void **)(&zoz));
29448a46eb9SPierre Jolivet       if (jac->constants[2]) PetscCallExternal(HYPRE_IJVectorGetObject, jac->constants[2]->ij, (void **)(&zzo));
295792fecdfSBarry Smith       PetscCallExternal(HYPRE_AMSSetEdgeConstantVectors, jac->hsolver, ozz, zoz, zzo);
2965ac14e1cSStefano Zampini     }
2975ac14e1cSStefano Zampini     if (jac->coords[0]) {
2985ac14e1cSStefano Zampini       HYPRE_ParVector coords[3];
2995ac14e1cSStefano Zampini       coords[0] = NULL;
3005ac14e1cSStefano Zampini       coords[1] = NULL;
3015ac14e1cSStefano Zampini       coords[2] = NULL;
302792fecdfSBarry Smith       if (jac->coords[0]) PetscCallExternal(HYPRE_IJVectorGetObject, jac->coords[0]->ij, (void **)(&coords[0]));
303792fecdfSBarry Smith       if (jac->coords[1]) PetscCallExternal(HYPRE_IJVectorGetObject, jac->coords[1]->ij, (void **)(&coords[1]));
304792fecdfSBarry Smith       if (jac->coords[2]) PetscCallExternal(HYPRE_IJVectorGetObject, jac->coords[2]->ij, (void **)(&coords[2]));
305792fecdfSBarry Smith       PetscCallExternal(HYPRE_AMSSetCoordinateVectors, jac->hsolver, coords[0], coords[1], coords[2]);
3065ac14e1cSStefano Zampini     }
3075f80ce2aSJacob Faibussowitsch     PetscCheck(jac->G, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "HYPRE AMS preconditioner needs the discrete gradient operator via PCHYPRESetDiscreteGradient");
3085ac14e1cSStefano Zampini     hm = (Mat_HYPRE *)(jac->G->data);
309792fecdfSBarry Smith     PetscCallExternal(HYPRE_IJMatrixGetObject, hm->ij, (void **)(&parcsr));
310792fecdfSBarry Smith     PetscCallExternal(HYPRE_AMSSetDiscreteGradient, jac->hsolver, parcsr);
3115ac14e1cSStefano Zampini     if (jac->alpha_Poisson) {
3125ac14e1cSStefano Zampini       hm = (Mat_HYPRE *)(jac->alpha_Poisson->data);
313792fecdfSBarry Smith       PetscCallExternal(HYPRE_IJMatrixGetObject, hm->ij, (void **)(&parcsr));
314792fecdfSBarry Smith       PetscCallExternal(HYPRE_AMSSetAlphaPoissonMatrix, jac->hsolver, parcsr);
3155ac14e1cSStefano Zampini     }
3165ac14e1cSStefano Zampini     if (jac->ams_beta_is_zero) {
317792fecdfSBarry Smith       PetscCallExternal(HYPRE_AMSSetBetaPoissonMatrix, jac->hsolver, NULL);
3185ac14e1cSStefano Zampini     } else if (jac->beta_Poisson) {
3195ac14e1cSStefano Zampini       hm = (Mat_HYPRE *)(jac->beta_Poisson->data);
320792fecdfSBarry Smith       PetscCallExternal(HYPRE_IJMatrixGetObject, hm->ij, (void **)(&parcsr));
321792fecdfSBarry Smith       PetscCallExternal(HYPRE_AMSSetBetaPoissonMatrix, jac->hsolver, parcsr);
322be14dc20SKerry Key     } else if (jac->ams_beta_is_zero_part) {
323be14dc20SKerry Key       if (jac->interior) {
324be14dc20SKerry Key         HYPRE_ParVector interior = NULL;
325be14dc20SKerry Key         PetscCallExternal(HYPRE_IJVectorGetObject, jac->interior->ij, (void **)(&interior));
326be14dc20SKerry Key         PetscCallExternal(HYPRE_AMSSetInteriorNodes, jac->hsolver, interior);
327be14dc20SKerry Key       } else {
328be14dc20SKerry Key         jac->ams_beta_is_zero_part = PETSC_FALSE;
329be14dc20SKerry Key       }
3305ac14e1cSStefano Zampini     }
3316bf688a0SCe Qin     if (jac->ND_PiFull || (jac->ND_Pi[0] && jac->ND_Pi[1])) {
3326bf688a0SCe Qin       PetscInt           i;
3336bf688a0SCe Qin       HYPRE_ParCSRMatrix nd_parcsrfull, nd_parcsr[3];
3346bf688a0SCe Qin       if (jac->ND_PiFull) {
3356bf688a0SCe Qin         hm = (Mat_HYPRE *)(jac->ND_PiFull->data);
336792fecdfSBarry Smith         PetscCallExternal(HYPRE_IJMatrixGetObject, hm->ij, (void **)(&nd_parcsrfull));
3376bf688a0SCe Qin       } else {
3386bf688a0SCe Qin         nd_parcsrfull = NULL;
3396bf688a0SCe Qin       }
3406bf688a0SCe Qin       for (i = 0; i < 3; ++i) {
3416bf688a0SCe Qin         if (jac->ND_Pi[i]) {
3426bf688a0SCe Qin           hm = (Mat_HYPRE *)(jac->ND_Pi[i]->data);
343792fecdfSBarry Smith           PetscCallExternal(HYPRE_IJMatrixGetObject, hm->ij, (void **)(&nd_parcsr[i]));
3446bf688a0SCe Qin         } else {
3456bf688a0SCe Qin           nd_parcsr[i] = NULL;
3466bf688a0SCe Qin         }
3476bf688a0SCe Qin       }
348792fecdfSBarry Smith       PetscCallExternal(HYPRE_AMSSetInterpolations, jac->hsolver, nd_parcsrfull, nd_parcsr[0], nd_parcsr[1], nd_parcsr[2]);
3496bf688a0SCe Qin     }
3504cb006feSStefano Zampini   }
351863406b8SStefano Zampini   /* special case for ADS */
352863406b8SStefano Zampini   if (jac->setup == HYPRE_ADSSetup) {
3535ac14e1cSStefano Zampini     Mat_HYPRE         *hm;
3545ac14e1cSStefano Zampini     HYPRE_ParCSRMatrix parcsr;
3556bf688a0SCe Qin     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])))) {
3566bf688a0SCe Qin       SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "HYPRE ADS preconditioner needs either the coordinate vectors via PCSetCoordinates() or the interpolation matrices via PCHYPRESetInterpolations");
3579371c9d4SSatish Balay     } 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");
3585f80ce2aSJacob Faibussowitsch     PetscCheck(jac->G, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "HYPRE ADS preconditioner needs the discrete gradient operator via PCHYPRESetDiscreteGradient");
3595f80ce2aSJacob Faibussowitsch     PetscCheck(jac->C, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "HYPRE ADS preconditioner needs the discrete curl operator via PCHYPRESetDiscreteGradient");
3605ac14e1cSStefano Zampini     if (jac->coords[0]) {
3615ac14e1cSStefano Zampini       HYPRE_ParVector coords[3];
3625ac14e1cSStefano Zampini       coords[0] = NULL;
3635ac14e1cSStefano Zampini       coords[1] = NULL;
3645ac14e1cSStefano Zampini       coords[2] = NULL;
365792fecdfSBarry Smith       if (jac->coords[0]) PetscCallExternal(HYPRE_IJVectorGetObject, jac->coords[0]->ij, (void **)(&coords[0]));
366792fecdfSBarry Smith       if (jac->coords[1]) PetscCallExternal(HYPRE_IJVectorGetObject, jac->coords[1]->ij, (void **)(&coords[1]));
367792fecdfSBarry Smith       if (jac->coords[2]) PetscCallExternal(HYPRE_IJVectorGetObject, jac->coords[2]->ij, (void **)(&coords[2]));
368792fecdfSBarry Smith       PetscCallExternal(HYPRE_ADSSetCoordinateVectors, jac->hsolver, coords[0], coords[1], coords[2]);
3695ac14e1cSStefano Zampini     }
3705ac14e1cSStefano Zampini     hm = (Mat_HYPRE *)(jac->G->data);
371792fecdfSBarry Smith     PetscCallExternal(HYPRE_IJMatrixGetObject, hm->ij, (void **)(&parcsr));
372792fecdfSBarry Smith     PetscCallExternal(HYPRE_ADSSetDiscreteGradient, jac->hsolver, parcsr);
3735ac14e1cSStefano Zampini     hm = (Mat_HYPRE *)(jac->C->data);
374792fecdfSBarry Smith     PetscCallExternal(HYPRE_IJMatrixGetObject, hm->ij, (void **)(&parcsr));
375792fecdfSBarry Smith     PetscCallExternal(HYPRE_ADSSetDiscreteCurl, jac->hsolver, parcsr);
3766bf688a0SCe Qin     if ((jac->RT_PiFull || (jac->RT_Pi[0] && jac->RT_Pi[1])) && (jac->ND_PiFull || (jac->ND_Pi[0] && jac->ND_Pi[1]))) {
3776bf688a0SCe Qin       PetscInt           i;
3786bf688a0SCe Qin       HYPRE_ParCSRMatrix rt_parcsrfull, rt_parcsr[3];
3796bf688a0SCe Qin       HYPRE_ParCSRMatrix nd_parcsrfull, nd_parcsr[3];
3806bf688a0SCe Qin       if (jac->RT_PiFull) {
3816bf688a0SCe Qin         hm = (Mat_HYPRE *)(jac->RT_PiFull->data);
382792fecdfSBarry Smith         PetscCallExternal(HYPRE_IJMatrixGetObject, hm->ij, (void **)(&rt_parcsrfull));
3836bf688a0SCe Qin       } else {
3846bf688a0SCe Qin         rt_parcsrfull = NULL;
3856bf688a0SCe Qin       }
3866bf688a0SCe Qin       for (i = 0; i < 3; ++i) {
3876bf688a0SCe Qin         if (jac->RT_Pi[i]) {
3886bf688a0SCe Qin           hm = (Mat_HYPRE *)(jac->RT_Pi[i]->data);
389792fecdfSBarry Smith           PetscCallExternal(HYPRE_IJMatrixGetObject, hm->ij, (void **)(&rt_parcsr[i]));
3906bf688a0SCe Qin         } else {
3916bf688a0SCe Qin           rt_parcsr[i] = NULL;
3926bf688a0SCe Qin         }
3936bf688a0SCe Qin       }
3946bf688a0SCe Qin       if (jac->ND_PiFull) {
3956bf688a0SCe Qin         hm = (Mat_HYPRE *)(jac->ND_PiFull->data);
396792fecdfSBarry Smith         PetscCallExternal(HYPRE_IJMatrixGetObject, hm->ij, (void **)(&nd_parcsrfull));
3976bf688a0SCe Qin       } else {
3986bf688a0SCe Qin         nd_parcsrfull = NULL;
3996bf688a0SCe Qin       }
4006bf688a0SCe Qin       for (i = 0; i < 3; ++i) {
4016bf688a0SCe Qin         if (jac->ND_Pi[i]) {
4026bf688a0SCe Qin           hm = (Mat_HYPRE *)(jac->ND_Pi[i]->data);
403792fecdfSBarry Smith           PetscCallExternal(HYPRE_IJMatrixGetObject, hm->ij, (void **)(&nd_parcsr[i]));
4046bf688a0SCe Qin         } else {
4056bf688a0SCe Qin           nd_parcsr[i] = NULL;
4066bf688a0SCe Qin         }
4076bf688a0SCe Qin       }
408792fecdfSBarry Smith       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]);
4096bf688a0SCe Qin     }
410863406b8SStefano Zampini   }
411792fecdfSBarry Smith   PetscCallExternal(HYPRE_IJMatrixGetObject, hjac->ij, (void **)&hmat);
412792fecdfSBarry Smith   PetscCallExternal(HYPRE_IJVectorGetObject, hjac->b->ij, (void **)&bv);
413792fecdfSBarry Smith   PetscCallExternal(HYPRE_IJVectorGetObject, hjac->x->ij, (void **)&xv);
414792fecdfSBarry Smith   PetscCallExternal(jac->setup, jac->hsolver, hmat, bv, xv);
41516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
41616d9e3a6SLisandro Dalcin }
41716d9e3a6SLisandro Dalcin 
4189371c9d4SSatish Balay static PetscErrorCode PCApply_HYPRE(PC pc, Vec b, Vec x) {
41916d9e3a6SLisandro Dalcin   PC_HYPRE          *jac  = (PC_HYPRE *)pc->data;
42049a781f5SStefano Zampini   Mat_HYPRE         *hjac = (Mat_HYPRE *)(jac->hpmat->data);
42116d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
42216d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv, jxv;
42316d9e3a6SLisandro Dalcin 
42416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
4259566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(hypreCitation, &cite));
4269566063dSJacob Faibussowitsch   if (!jac->applyrichardson) PetscCall(VecSet(x, 0.0));
4279566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorPushVecRead(hjac->b, b));
4289566063dSJacob Faibussowitsch   if (jac->applyrichardson) PetscCall(VecHYPRE_IJVectorPushVec(hjac->x, x));
4299566063dSJacob Faibussowitsch   else PetscCall(VecHYPRE_IJVectorPushVecWrite(hjac->x, x));
430792fecdfSBarry Smith   PetscCallExternal(HYPRE_IJMatrixGetObject, hjac->ij, (void **)&hmat);
431792fecdfSBarry Smith   PetscCallExternal(HYPRE_IJVectorGetObject, hjac->b->ij, (void **)&jbv);
432792fecdfSBarry Smith   PetscCallExternal(HYPRE_IJVectorGetObject, hjac->x->ij, (void **)&jxv);
4339371c9d4SSatish Balay   PetscStackCallExternalVoid(
4349371c9d4SSatish Balay     "Hypre solve", do {
4355f80ce2aSJacob Faibussowitsch       HYPRE_Int hierr = (*jac->solve)(jac->hsolver, hmat, jbv, jxv);
4365f80ce2aSJacob Faibussowitsch       if (hierr) {
4375f80ce2aSJacob Faibussowitsch         PetscCheck(hierr == HYPRE_ERROR_CONV, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in HYPRE solver, error code %d", (int)hierr);
4385f80ce2aSJacob Faibussowitsch         hypre__global_error = 0;
4395f80ce2aSJacob Faibussowitsch       }
4405f80ce2aSJacob Faibussowitsch     } while (0));
44116d9e3a6SLisandro Dalcin 
44248a46eb9SPierre Jolivet   if (jac->setup == HYPRE_AMSSetup && jac->ams_beta_is_zero_part) PetscCallExternal(HYPRE_AMSProjectOutGradients, jac->hsolver, jxv);
4439566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorPopVec(hjac->x));
4449566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorPopVec(hjac->b));
44516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
44616d9e3a6SLisandro Dalcin }
44716d9e3a6SLisandro Dalcin 
4489371c9d4SSatish Balay static PetscErrorCode PCReset_HYPRE(PC pc) {
4498695de01SBarry Smith   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
4508695de01SBarry Smith 
4518695de01SBarry Smith   PetscFunctionBegin;
4529566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->hpmat));
4539566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->G));
4549566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->C));
4559566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->alpha_Poisson));
4569566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->beta_Poisson));
4579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->RT_PiFull));
4589566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->RT_Pi[0]));
4599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->RT_Pi[1]));
4609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->RT_Pi[2]));
4619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->ND_PiFull));
4629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->ND_Pi[0]));
4639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->ND_Pi[1]));
4649566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->ND_Pi[2]));
4659566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorDestroy(&jac->coords[0]));
4669566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorDestroy(&jac->coords[1]));
4679566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorDestroy(&jac->coords[2]));
4689566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorDestroy(&jac->constants[0]));
4699566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorDestroy(&jac->constants[1]));
4709566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorDestroy(&jac->constants[2]));
471be14dc20SKerry Key   PetscCall(VecHYPRE_IJVectorDestroy(&jac->interior));
4729566063dSJacob Faibussowitsch   PetscCall(PCHYPREResetNearNullSpace_Private(pc));
4735ac14e1cSStefano Zampini   jac->ams_beta_is_zero      = PETSC_FALSE;
474be14dc20SKerry Key   jac->ams_beta_is_zero_part = PETSC_FALSE;
4755ac14e1cSStefano Zampini   jac->dim                   = 0;
4768695de01SBarry Smith   PetscFunctionReturn(0);
4778695de01SBarry Smith }
4788695de01SBarry Smith 
4799371c9d4SSatish Balay static PetscErrorCode PCDestroy_HYPRE(PC pc) {
48016d9e3a6SLisandro Dalcin   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
48116d9e3a6SLisandro Dalcin 
48216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
4839566063dSJacob Faibussowitsch   PetscCall(PCReset_HYPRE(pc));
484792fecdfSBarry Smith   if (jac->destroy) PetscCallExternal(jac->destroy, jac->hsolver);
4859566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->hypre_type));
486db6f9c32SMark Adams #if PETSC_PKG_HYPRE_VERSION_GE(2, 23, 0)
4879566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->spgemm_type));
488db6f9c32SMark Adams #endif
4899566063dSJacob Faibussowitsch   if (jac->comm_hypre != MPI_COMM_NULL) PetscCall(PetscCommRestoreComm(PetscObjectComm((PetscObject)pc), &jac->comm_hypre));
4909566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
49116d9e3a6SLisandro Dalcin 
4929566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, 0));
4939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPRESetType_C", NULL));
4949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPREGetType_C", NULL));
4959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPRESetDiscreteGradient_C", NULL));
4969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPRESetDiscreteCurl_C", NULL));
4979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPRESetInterpolations_C", NULL));
4989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPRESetConstantEdgeVectors_C", NULL));
4999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPRESetPoissonMatrix_C", NULL));
5002e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPRESetEdgeConstantVectors_C", NULL));
501be14dc20SKerry Key   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPREAMSSetInteriorNodes_C", NULL));
5029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", NULL));
5039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", NULL));
5049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGalerkinSetMatProductAlgorithm_C", NULL));
5059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGalerkinGetMatProductAlgorithm_C", NULL));
5062e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", NULL));
50716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
50816d9e3a6SLisandro Dalcin }
50916d9e3a6SLisandro Dalcin 
5109371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc, PetscOptionItems *PetscOptionsObject) {
51116d9e3a6SLisandro Dalcin   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
512ace3abfcSBarry Smith   PetscBool flag;
51316d9e3a6SLisandro Dalcin 
51416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
515d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "HYPRE Pilut Options");
5169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_pilut_maxiter", "Number of iterations", "None", jac->maxiter, &jac->maxiter, &flag));
517792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_ParCSRPilutSetMaxIter, jac->hsolver, jac->maxiter);
5189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_pilut_tol", "Drop tolerance", "None", jac->tol, &jac->tol, &flag));
519792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_ParCSRPilutSetDropTolerance, jac->hsolver, jac->tol);
5209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_pilut_factorrowsize", "FactorRowSize", "None", jac->factorrowsize, &jac->factorrowsize, &flag));
521792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_ParCSRPilutSetFactorRowSize, jac->hsolver, jac->factorrowsize);
522d0609cedSBarry Smith   PetscOptionsHeadEnd();
52316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
52416d9e3a6SLisandro Dalcin }
52516d9e3a6SLisandro Dalcin 
5269371c9d4SSatish Balay static PetscErrorCode PCView_HYPRE_Pilut(PC pc, PetscViewer viewer) {
52716d9e3a6SLisandro Dalcin   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
528ace3abfcSBarry Smith   PetscBool iascii;
52916d9e3a6SLisandro Dalcin 
53016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
5319566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
53216d9e3a6SLisandro Dalcin   if (iascii) {
5339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  HYPRE Pilut preconditioning\n"));
53416d9e3a6SLisandro Dalcin     if (jac->maxiter != PETSC_DEFAULT) {
53563a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    maximum number of iterations %" PetscInt_FMT "\n", jac->maxiter));
53616d9e3a6SLisandro Dalcin     } else {
5379566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    default maximum number of iterations \n"));
53816d9e3a6SLisandro Dalcin     }
53916d9e3a6SLisandro Dalcin     if (jac->tol != PETSC_DEFAULT) {
5409566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    drop tolerance %g\n", (double)jac->tol));
54116d9e3a6SLisandro Dalcin     } else {
5429566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    default drop tolerance \n"));
54316d9e3a6SLisandro Dalcin     }
54416d9e3a6SLisandro Dalcin     if (jac->factorrowsize != PETSC_DEFAULT) {
54563a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    factor row size %" PetscInt_FMT "\n", jac->factorrowsize));
54616d9e3a6SLisandro Dalcin     } else {
5479566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    default factor row size \n"));
54816d9e3a6SLisandro Dalcin     }
54916d9e3a6SLisandro Dalcin   }
55016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
55116d9e3a6SLisandro Dalcin }
55216d9e3a6SLisandro Dalcin 
5539371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc, PetscOptionItems *PetscOptionsObject) {
554db966c6cSHong Zhang   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
5558bf83915SBarry Smith   PetscBool flag, eu_bj = jac->eu_bj ? PETSC_TRUE : PETSC_FALSE;
556db966c6cSHong Zhang 
557db966c6cSHong Zhang   PetscFunctionBegin;
558d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "HYPRE Euclid Options");
5599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_euclid_level", "Factorization levels", "None", jac->eu_level, &jac->eu_level, &flag));
560792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_EuclidSetLevel, jac->hsolver, jac->eu_level);
5618bf83915SBarry Smith 
5629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_euclid_droptolerance", "Drop tolerance for ILU(k) in Euclid", "None", jac->eu_droptolerance, &jac->eu_droptolerance, &flag));
5638bf83915SBarry Smith   if (flag) {
5648bf83915SBarry Smith     PetscMPIInt size;
5658bf83915SBarry Smith 
5669566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
5677827d75bSBarry Smith     PetscCheck(size == 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "hypre's Euclid does not support a parallel drop tolerance");
568792fecdfSBarry Smith     PetscCallExternal(HYPRE_EuclidSetILUT, jac->hsolver, jac->eu_droptolerance);
5698bf83915SBarry Smith   }
5708bf83915SBarry Smith 
5719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_hypre_euclid_bj", "Use Block Jacobi for ILU in Euclid", "None", eu_bj, &eu_bj, &flag));
5728bf83915SBarry Smith   if (flag) {
5738bf83915SBarry Smith     jac->eu_bj = eu_bj ? 1 : 0;
574792fecdfSBarry Smith     PetscCallExternal(HYPRE_EuclidSetBJ, jac->hsolver, jac->eu_bj);
5758bf83915SBarry Smith   }
576d0609cedSBarry Smith   PetscOptionsHeadEnd();
577db966c6cSHong Zhang   PetscFunctionReturn(0);
578db966c6cSHong Zhang }
579db966c6cSHong Zhang 
5809371c9d4SSatish Balay static PetscErrorCode PCView_HYPRE_Euclid(PC pc, PetscViewer viewer) {
581db966c6cSHong Zhang   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
582db966c6cSHong Zhang   PetscBool iascii;
583db966c6cSHong Zhang 
584db966c6cSHong Zhang   PetscFunctionBegin;
5859566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
586db966c6cSHong Zhang   if (iascii) {
5879566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  HYPRE Euclid preconditioning\n"));
588db966c6cSHong Zhang     if (jac->eu_level != PETSC_DEFAULT) {
58963a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    factorization levels %" PetscInt_FMT "\n", jac->eu_level));
590db966c6cSHong Zhang     } else {
5919566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    default factorization levels \n"));
592db966c6cSHong Zhang     }
5939566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    drop tolerance %g\n", (double)jac->eu_droptolerance));
59463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    use Block-Jacobi? %" PetscInt_FMT "\n", jac->eu_bj));
595db966c6cSHong Zhang   }
596db966c6cSHong Zhang   PetscFunctionReturn(0);
597db966c6cSHong Zhang }
598db966c6cSHong Zhang 
5999371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc, Vec b, Vec x) {
60016d9e3a6SLisandro Dalcin   PC_HYPRE          *jac  = (PC_HYPRE *)pc->data;
60149a781f5SStefano Zampini   Mat_HYPRE         *hjac = (Mat_HYPRE *)(jac->hpmat->data);
60216d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
60316d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv, jxv;
60416d9e3a6SLisandro Dalcin 
60516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
6069566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(hypreCitation, &cite));
6079566063dSJacob Faibussowitsch   PetscCall(VecSet(x, 0.0));
6089566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorPushVecRead(hjac->x, b));
6099566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorPushVecWrite(hjac->b, x));
61016d9e3a6SLisandro Dalcin 
611792fecdfSBarry Smith   PetscCallExternal(HYPRE_IJMatrixGetObject, hjac->ij, (void **)&hmat);
612792fecdfSBarry Smith   PetscCallExternal(HYPRE_IJVectorGetObject, hjac->b->ij, (void **)&jbv);
613792fecdfSBarry Smith   PetscCallExternal(HYPRE_IJVectorGetObject, hjac->x->ij, (void **)&jxv);
61416d9e3a6SLisandro Dalcin 
6159371c9d4SSatish Balay   PetscStackCallExternalVoid(
6169371c9d4SSatish Balay     "Hypre Transpose solve", do {
6175f80ce2aSJacob Faibussowitsch       HYPRE_Int hierr = HYPRE_BoomerAMGSolveT(jac->hsolver, hmat, jbv, jxv);
6185f80ce2aSJacob Faibussowitsch       if (hierr) {
61916d9e3a6SLisandro Dalcin         /* error code of 1 in BoomerAMG merely means convergence not achieved */
6205f80ce2aSJacob Faibussowitsch         PetscCheck(hierr == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in HYPRE solver, error code %d", (int)hierr);
6215f80ce2aSJacob Faibussowitsch         hypre__global_error = 0;
6225f80ce2aSJacob Faibussowitsch       }
6235f80ce2aSJacob Faibussowitsch     } while (0));
62416d9e3a6SLisandro Dalcin 
6259566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorPopVec(hjac->x));
6269566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorPopVec(hjac->b));
62716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
62816d9e3a6SLisandro Dalcin }
62916d9e3a6SLisandro Dalcin 
6309371c9d4SSatish Balay static PetscErrorCode PCMGGalerkinSetMatProductAlgorithm_HYPRE_BoomerAMG(PC pc, const char name[]) {
631db6f9c32SMark Adams   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
632db6f9c32SMark Adams   PetscBool flag;
633db6f9c32SMark Adams 
634db6f9c32SMark Adams #if PETSC_PKG_HYPRE_VERSION_GE(2, 23, 0)
635db6f9c32SMark Adams   PetscFunctionBegin;
636db6f9c32SMark Adams   if (jac->spgemm_type) {
6379566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(jac->spgemm_type, name, &flag));
63828b400f6SJacob Faibussowitsch     PetscCheck(flag, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Cannot reset the HYPRE SpGEMM (really we can)");
639db6f9c32SMark Adams     PetscFunctionReturn(0);
640db6f9c32SMark Adams   } else {
6419566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(name, &jac->spgemm_type));
642db6f9c32SMark Adams   }
6439566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("cusparse", jac->spgemm_type, &flag));
644db6f9c32SMark Adams   if (flag) {
645792fecdfSBarry Smith     PetscCallExternal(HYPRE_SetSpGemmUseCusparse, 1);
646db6f9c32SMark Adams     PetscFunctionReturn(0);
647db6f9c32SMark Adams   }
6489566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("hypre", jac->spgemm_type, &flag));
649db6f9c32SMark Adams   if (flag) {
650792fecdfSBarry Smith     PetscCallExternal(HYPRE_SetSpGemmUseCusparse, 0);
651db6f9c32SMark Adams     PetscFunctionReturn(0);
652db6f9c32SMark Adams   }
653db6f9c32SMark Adams   jac->spgemm_type = NULL;
65498921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown HYPRE SpGEM type %s; Choices are cusparse, hypre", name);
655db6f9c32SMark Adams #endif
656db6f9c32SMark Adams }
657db6f9c32SMark Adams 
6589371c9d4SSatish Balay static PetscErrorCode PCMGGalerkinGetMatProductAlgorithm_HYPRE_BoomerAMG(PC pc, const char *spgemm[]) {
659db6f9c32SMark Adams   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
660db6f9c32SMark Adams 
661db6f9c32SMark Adams   PetscFunctionBegin;
662db6f9c32SMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
663db6f9c32SMark Adams #if PETSC_PKG_HYPRE_VERSION_GE(2, 23, 0)
664db6f9c32SMark Adams   *spgemm = jac->spgemm_type;
665db6f9c32SMark Adams #endif
666db6f9c32SMark Adams   PetscFunctionReturn(0);
667db6f9c32SMark Adams }
668db6f9c32SMark Adams 
66916d9e3a6SLisandro Dalcin static const char    *HYPREBoomerAMGCycleType[]   = {"", "V", "W"};
6700f1074feSSatish Balay static const char    *HYPREBoomerAMGCoarsenType[] = {"CLJP", "Ruge-Stueben", "", "modifiedRuge-Stueben", "", "", "Falgout", "", "PMIS", "", "HMIS"};
67116d9e3a6SLisandro Dalcin static const char    *HYPREBoomerAMGMeasureType[] = {"local", "global"};
67265de4495SJed Brown /* The following corresponds to HYPRE_BoomerAMGSetRelaxType which has many missing numbers in the enum */
6736a251517SEike Mueller static const char    *HYPREBoomerAMGSmoothType[]  = {"Schwarz-smoothers", "Pilut", "ParaSails", "Euclid"};
6749371c9d4SSatish Balay static const char    *HYPREBoomerAMGRelaxType[]   = {"Jacobi", "sequential-Gauss-Seidel", "seqboundary-Gauss-Seidel", "SOR/Jacobi", "backward-SOR/Jacobi", "" /* [5] hybrid chaotic Gauss-Seidel (works only with OpenMP) */, "symmetric-SOR/Jacobi", "" /* 7 */, "l1scaled-SOR/Jacobi", "Gaussian-elimination", "" /* 10 */, "" /* 11 */, "" /* 12 */, "l1-Gauss-Seidel" /* nonsymmetric */, "backward-l1-Gauss-Seidel" /* nonsymmetric */, "CG" /* non-stationary */, "Chebyshev", "FCF-Jacobi", "l1scaled-Jacobi"};
6759371c9d4SSatish Balay static const char    *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i", "ext+i-cc", "standard", "standard-wts", "block", "block-wtd", "FF", "FF1", "ext", "ad-wts", "ext-mm", "ext+i-mm", "ext+e-mm"};
6769371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc, PetscOptionItems *PetscOptionsObject) {
67716d9e3a6SLisandro Dalcin   PC_HYPRE   *jac = (PC_HYPRE *)pc->data;
67822e51d31SStefano Zampini   PetscInt    bs, n, indx, level;
679ace3abfcSBarry Smith   PetscBool   flg, tmp_truth;
68016d9e3a6SLisandro Dalcin   double      tmpdbl, twodbl[2];
681589dcaf0SStefano Zampini   const char *symtlist[]           = {"nonsymmetric", "SPD", "nonsymmetric,SPD"};
682db6f9c32SMark Adams   const char *PCHYPRESpgemmTypes[] = {"cusparse", "hypre"};
68316d9e3a6SLisandro Dalcin 
68416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
685d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "HYPRE BoomerAMG Options");
6869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_cycle_type", "Cycle type", "None", HYPREBoomerAMGCycleType + 1, 2, HYPREBoomerAMGCycleType[jac->cycletype], &indx, &flg));
68716d9e3a6SLisandro Dalcin   if (flg) {
6884336a9eeSBarry Smith     jac->cycletype = indx + 1;
689792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetCycleType, jac->hsolver, jac->cycletype);
69016d9e3a6SLisandro Dalcin   }
6919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_max_levels", "Number of levels (of grids) allowed", "None", jac->maxlevels, &jac->maxlevels, &flg));
69216d9e3a6SLisandro Dalcin   if (flg) {
69363a3b9bcSJacob Faibussowitsch     PetscCheck(jac->maxlevels >= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Number of levels %" PetscInt_FMT " must be at least two", jac->maxlevels);
694792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetMaxLevels, jac->hsolver, jac->maxlevels);
69516d9e3a6SLisandro Dalcin   }
6969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_max_iter", "Maximum iterations used PER hypre call", "None", jac->maxiter, &jac->maxiter, &flg));
69716d9e3a6SLisandro Dalcin   if (flg) {
69863a3b9bcSJacob Faibussowitsch     PetscCheck(jac->maxiter >= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Number of iterations %" PetscInt_FMT " must be at least one", jac->maxiter);
699792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetMaxIter, jac->hsolver, jac->maxiter);
70016d9e3a6SLisandro Dalcin   }
7019566063dSJacob Faibussowitsch   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));
70216d9e3a6SLisandro Dalcin   if (flg) {
70308401ef6SPierre Jolivet     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);
704792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetTol, jac->hsolver, jac->tol);
70516d9e3a6SLisandro Dalcin   }
70622e51d31SStefano Zampini   bs = 1;
70748a46eb9SPierre Jolivet   if (pc->pmat) PetscCall(MatGetBlockSize(pc->pmat, &bs));
7089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_numfunctions", "Number of functions", "HYPRE_BoomerAMGSetNumFunctions", bs, &bs, &flg));
70948a46eb9SPierre Jolivet   if (flg) PetscCallExternal(HYPRE_BoomerAMGSetNumFunctions, jac->hsolver, bs);
71016d9e3a6SLisandro Dalcin 
7119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_boomeramg_truncfactor", "Truncation factor for interpolation (0=no truncation)", "None", jac->truncfactor, &jac->truncfactor, &flg));
71216d9e3a6SLisandro Dalcin   if (flg) {
71308401ef6SPierre Jolivet     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);
714792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetTruncFactor, jac->hsolver, jac->truncfactor);
71516d9e3a6SLisandro Dalcin   }
71616d9e3a6SLisandro Dalcin 
7179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_P_max", "Max elements per row for interpolation operator (0=unlimited)", "None", jac->pmax, &jac->pmax, &flg));
7180f1074feSSatish Balay   if (flg) {
71963a3b9bcSJacob Faibussowitsch     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);
720792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetPMaxElmts, jac->hsolver, jac->pmax);
7210f1074feSSatish Balay   }
7220f1074feSSatish Balay 
7239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-pc_hypre_boomeramg_agg_nl", "Number of levels of aggressive coarsening", "None", jac->agg_nl, &jac->agg_nl, &flg, 0, jac->maxlevels));
724792fecdfSBarry Smith   if (flg) PetscCallExternal(HYPRE_BoomerAMGSetAggNumLevels, jac->hsolver, jac->agg_nl);
7250f1074feSSatish Balay 
7269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_agg_num_paths", "Number of paths for aggressive coarsening", "None", jac->agg_num_paths, &jac->agg_num_paths, &flg));
7270f1074feSSatish Balay   if (flg) {
72863a3b9bcSJacob Faibussowitsch     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);
729792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetNumPaths, jac->hsolver, jac->agg_num_paths);
7300f1074feSSatish Balay   }
7310f1074feSSatish Balay 
7329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_boomeramg_strong_threshold", "Threshold for being strongly connected", "None", jac->strongthreshold, &jac->strongthreshold, &flg));
73316d9e3a6SLisandro Dalcin   if (flg) {
73408401ef6SPierre Jolivet     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);
735792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetStrongThreshold, jac->hsolver, jac->strongthreshold);
73616d9e3a6SLisandro Dalcin   }
7379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_boomeramg_max_row_sum", "Maximum row sum", "None", jac->maxrowsum, &jac->maxrowsum, &flg));
73816d9e3a6SLisandro Dalcin   if (flg) {
73908401ef6SPierre Jolivet     PetscCheck(jac->maxrowsum >= 0.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Maximum row sum %g must be greater than zero", (double)jac->maxrowsum);
74008401ef6SPierre Jolivet     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);
741792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetMaxRowSum, jac->hsolver, jac->maxrowsum);
74216d9e3a6SLisandro Dalcin   }
74316d9e3a6SLisandro Dalcin 
74416d9e3a6SLisandro Dalcin   /* Grid sweeps */
7459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_all", "Number of sweeps for the up and down grid levels", "None", jac->gridsweeps[0], &indx, &flg));
74616d9e3a6SLisandro Dalcin   if (flg) {
747792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetNumSweeps, jac->hsolver, indx);
74816d9e3a6SLisandro Dalcin     /* modify the jac structure so we can view the updated options with PC_View */
74916d9e3a6SLisandro Dalcin     jac->gridsweeps[0] = indx;
7500f1074feSSatish Balay     jac->gridsweeps[1] = indx;
7510f1074feSSatish Balay     /*defaults coarse to 1 */
7520f1074feSSatish Balay     jac->gridsweeps[2] = 1;
75316d9e3a6SLisandro Dalcin   }
7549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_nodal_coarsen", "Use a nodal based coarsening 1-6", "HYPRE_BoomerAMGSetNodal", jac->nodal_coarsening, &jac->nodal_coarsening, &flg));
75548a46eb9SPierre Jolivet   if (flg) PetscCallExternal(HYPRE_BoomerAMGSetNodal, jac->hsolver, jac->nodal_coarsening);
7569566063dSJacob Faibussowitsch   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));
75748a46eb9SPierre Jolivet   if (flg) PetscCallExternal(HYPRE_BoomerAMGSetNodalDiag, jac->hsolver, jac->nodal_coarsening_diag);
7589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_vec_interp_variant", "Variant of algorithm 1-3", "HYPRE_BoomerAMGSetInterpVecVariant", jac->vec_interp_variant, &jac->vec_interp_variant, &flg));
75948a46eb9SPierre Jolivet   if (flg) PetscCallExternal(HYPRE_BoomerAMGSetInterpVecVariant, jac->hsolver, jac->vec_interp_variant);
7609566063dSJacob Faibussowitsch   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));
76148a46eb9SPierre Jolivet   if (flg) PetscCallExternal(HYPRE_BoomerAMGSetInterpVecQMax, jac->hsolver, jac->vec_interp_qmax);
7629566063dSJacob Faibussowitsch   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));
76348a46eb9SPierre Jolivet   if (flg) PetscCallExternal(HYPRE_BoomerAMGSetSmoothInterpVectors, jac->hsolver, jac->vec_interp_smooth);
7649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_interp_refine", "Preprocess the interpolation matrix through iterative weight refinement", "HYPRE_BoomerAMGSetInterpRefine", jac->interp_refine, &jac->interp_refine, &flg));
76548a46eb9SPierre Jolivet   if (flg) PetscCallExternal(HYPRE_BoomerAMGSetInterpRefine, jac->hsolver, jac->interp_refine);
7669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down", "Number of sweeps for the down cycles", "None", jac->gridsweeps[0], &indx, &flg));
76716d9e3a6SLisandro Dalcin   if (flg) {
768792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetCycleNumSweeps, jac->hsolver, indx, 1);
7690f1074feSSatish Balay     jac->gridsweeps[0] = indx;
77016d9e3a6SLisandro Dalcin   }
7719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up", "Number of sweeps for the up cycles", "None", jac->gridsweeps[1], &indx, &flg));
77216d9e3a6SLisandro Dalcin   if (flg) {
773792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetCycleNumSweeps, jac->hsolver, indx, 2);
7740f1074feSSatish Balay     jac->gridsweeps[1] = indx;
77516d9e3a6SLisandro Dalcin   }
7769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse", "Number of sweeps for the coarse level", "None", jac->gridsweeps[2], &indx, &flg));
77716d9e3a6SLisandro Dalcin   if (flg) {
778792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetCycleNumSweeps, jac->hsolver, indx, 3);
7790f1074feSSatish Balay     jac->gridsweeps[2] = indx;
78016d9e3a6SLisandro Dalcin   }
78116d9e3a6SLisandro Dalcin 
7826a251517SEike Mueller   /* Smooth type */
783dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_smooth_type", "Enable more complex smoothers", "None", HYPREBoomerAMGSmoothType, PETSC_STATIC_ARRAY_LENGTH(HYPREBoomerAMGSmoothType), HYPREBoomerAMGSmoothType[0], &indx, &flg));
7846a251517SEike Mueller   if (flg) {
7856a251517SEike Mueller     jac->smoothtype = indx;
786792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetSmoothType, jac->hsolver, indx + 6);
7878131ecf7SEike Mueller     jac->smoothnumlevels = 25;
788792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetSmoothNumLevels, jac->hsolver, 25);
7898131ecf7SEike Mueller   }
7908131ecf7SEike Mueller 
7918131ecf7SEike Mueller   /* Number of smoothing levels */
7929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_smooth_num_levels", "Number of levels on which more complex smoothers are used", "None", 25, &indx, &flg));
7938131ecf7SEike Mueller   if (flg && (jac->smoothtype != -1)) {
7948131ecf7SEike Mueller     jac->smoothnumlevels = indx;
795792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetSmoothNumLevels, jac->hsolver, indx);
7966a251517SEike Mueller   }
7976a251517SEike Mueller 
7981810e44eSEike Mueller   /* Number of levels for ILU(k) for Euclid */
7999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_eu_level", "Number of levels for ILU(k) in Euclid smoother", "None", 0, &indx, &flg));
8001810e44eSEike Mueller   if (flg && (jac->smoothtype == 3)) {
8011810e44eSEike Mueller     jac->eu_level = indx;
802792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetEuLevel, jac->hsolver, indx);
8031810e44eSEike Mueller   }
8041810e44eSEike Mueller 
8051810e44eSEike Mueller   /* Filter for ILU(k) for Euclid */
8061810e44eSEike Mueller   double droptolerance;
8079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_boomeramg_eu_droptolerance", "Drop tolerance for ILU(k) in Euclid smoother", "None", 0, &droptolerance, &flg));
8081810e44eSEike Mueller   if (flg && (jac->smoothtype == 3)) {
8091810e44eSEike Mueller     jac->eu_droptolerance = droptolerance;
810792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetEuLevel, jac->hsolver, droptolerance);
8111810e44eSEike Mueller   }
8121810e44eSEike Mueller 
8131810e44eSEike Mueller   /* Use Block Jacobi ILUT for Euclid */
8149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_hypre_boomeramg_eu_bj", "Use Block Jacobi for ILU in Euclid smoother?", "None", PETSC_FALSE, &tmp_truth, &flg));
8151810e44eSEike Mueller   if (flg && (jac->smoothtype == 3)) {
8161810e44eSEike Mueller     jac->eu_bj = tmp_truth;
817792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetEuBJ, jac->hsolver, jac->eu_bj);
8181810e44eSEike Mueller   }
8191810e44eSEike Mueller 
82016d9e3a6SLisandro Dalcin   /* Relax type */
821dd39110bSPierre Jolivet   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));
82216d9e3a6SLisandro Dalcin   if (flg) {
8230f1074feSSatish Balay     jac->relaxtype[0] = jac->relaxtype[1] = indx;
824792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetRelaxType, jac->hsolver, indx);
8250f1074feSSatish Balay     /* by default, coarse type set to 9 */
8260f1074feSSatish Balay     jac->relaxtype[2] = 9;
827792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetCycleRelaxType, jac->hsolver, 9, 3);
82816d9e3a6SLisandro Dalcin   }
829dd39110bSPierre Jolivet   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));
83016d9e3a6SLisandro Dalcin   if (flg) {
83116d9e3a6SLisandro Dalcin     jac->relaxtype[0] = indx;
832792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetCycleRelaxType, jac->hsolver, indx, 1);
83316d9e3a6SLisandro Dalcin   }
834dd39110bSPierre Jolivet   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));
83516d9e3a6SLisandro Dalcin   if (flg) {
8360f1074feSSatish Balay     jac->relaxtype[1] = indx;
837792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetCycleRelaxType, jac->hsolver, indx, 2);
83816d9e3a6SLisandro Dalcin   }
839dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse", "Relax type on coarse grid", "None", HYPREBoomerAMGRelaxType, PETSC_STATIC_ARRAY_LENGTH(HYPREBoomerAMGRelaxType), HYPREBoomerAMGRelaxType[9], &indx, &flg));
84016d9e3a6SLisandro Dalcin   if (flg) {
8410f1074feSSatish Balay     jac->relaxtype[2] = indx;
842792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetCycleRelaxType, jac->hsolver, indx, 3);
84316d9e3a6SLisandro Dalcin   }
84416d9e3a6SLisandro Dalcin 
84516d9e3a6SLisandro Dalcin   /* Relaxation Weight */
8469566063dSJacob Faibussowitsch   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));
84716d9e3a6SLisandro Dalcin   if (flg) {
848792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetRelaxWt, jac->hsolver, tmpdbl);
84916d9e3a6SLisandro Dalcin     jac->relaxweight = tmpdbl;
85016d9e3a6SLisandro Dalcin   }
85116d9e3a6SLisandro Dalcin 
85216d9e3a6SLisandro Dalcin   n         = 2;
85316d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
8549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-pc_hypre_boomeramg_relax_weight_level", "Set the relaxation weight for a particular level (weight,level)", "None", twodbl, &n, &flg));
85516d9e3a6SLisandro Dalcin   if (flg) {
85616d9e3a6SLisandro Dalcin     if (n == 2) {
85716d9e3a6SLisandro Dalcin       indx = (int)PetscAbsReal(twodbl[1]);
858792fecdfSBarry Smith       PetscCallExternal(HYPRE_BoomerAMGSetLevelRelaxWt, jac->hsolver, twodbl[0], indx);
85963a3b9bcSJacob Faibussowitsch     } 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);
86016d9e3a6SLisandro Dalcin   }
86116d9e3a6SLisandro Dalcin 
86216d9e3a6SLisandro Dalcin   /* Outer relaxation Weight */
8639566063dSJacob Faibussowitsch   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));
86416d9e3a6SLisandro Dalcin   if (flg) {
865792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetOuterWt, jac->hsolver, tmpdbl);
86616d9e3a6SLisandro Dalcin     jac->outerrelaxweight = tmpdbl;
86716d9e3a6SLisandro Dalcin   }
86816d9e3a6SLisandro Dalcin 
86916d9e3a6SLisandro Dalcin   n         = 2;
87016d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
8719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-pc_hypre_boomeramg_outer_relax_weight_level", "Set the outer relaxation weight for a particular level (weight,level)", "None", twodbl, &n, &flg));
87216d9e3a6SLisandro Dalcin   if (flg) {
87316d9e3a6SLisandro Dalcin     if (n == 2) {
87416d9e3a6SLisandro Dalcin       indx = (int)PetscAbsReal(twodbl[1]);
875792fecdfSBarry Smith       PetscCallExternal(HYPRE_BoomerAMGSetLevelOuterWt, jac->hsolver, twodbl[0], indx);
87663a3b9bcSJacob Faibussowitsch     } 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);
87716d9e3a6SLisandro Dalcin   }
87816d9e3a6SLisandro Dalcin 
87916d9e3a6SLisandro Dalcin   /* the Relax Order */
8809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg));
88116d9e3a6SLisandro Dalcin 
8828afaa268SBarry Smith   if (flg && tmp_truth) {
88316d9e3a6SLisandro Dalcin     jac->relaxorder = 0;
884792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetRelaxOrder, jac->hsolver, jac->relaxorder);
88516d9e3a6SLisandro Dalcin   }
886dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_measure_type", "Measure type", "None", HYPREBoomerAMGMeasureType, PETSC_STATIC_ARRAY_LENGTH(HYPREBoomerAMGMeasureType), HYPREBoomerAMGMeasureType[0], &indx, &flg));
88716d9e3a6SLisandro Dalcin   if (flg) {
88816d9e3a6SLisandro Dalcin     jac->measuretype = indx;
889792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetMeasureType, jac->hsolver, jac->measuretype);
89016d9e3a6SLisandro Dalcin   }
8910f1074feSSatish Balay   /* update list length 3/07 */
892dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type", "Coarsen type", "None", HYPREBoomerAMGCoarsenType, PETSC_STATIC_ARRAY_LENGTH(HYPREBoomerAMGCoarsenType), HYPREBoomerAMGCoarsenType[6], &indx, &flg));
89316d9e3a6SLisandro Dalcin   if (flg) {
89416d9e3a6SLisandro Dalcin     jac->coarsentype = indx;
895792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetCoarsenType, jac->hsolver, jac->coarsentype);
89616d9e3a6SLisandro Dalcin   }
8970f1074feSSatish Balay 
8989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_max_coarse_size", "Maximum size of coarsest grid", "None", jac->maxc, &jac->maxc, &flg));
89948a46eb9SPierre Jolivet   if (flg) PetscCallExternal(HYPRE_BoomerAMGSetMaxCoarseSize, jac->hsolver, jac->maxc);
9009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_min_coarse_size", "Minimum size of coarsest grid", "None", jac->minc, &jac->minc, &flg));
90148a46eb9SPierre Jolivet   if (flg) PetscCallExternal(HYPRE_BoomerAMGSetMinCoarseSize, jac->hsolver, jac->minc);
902db6f9c32SMark Adams #if PETSC_PKG_HYPRE_VERSION_GE(2, 23, 0)
903db6f9c32SMark Adams   // global parameter but is closely associated with BoomerAMG
904dd39110bSPierre Jolivet   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));
905db6f9c32SMark Adams   if (!flg) indx = 0;
9069566063dSJacob Faibussowitsch   PetscCall(PCMGGalerkinSetMatProductAlgorithm_HYPRE_BoomerAMG(pc, PCHYPRESpgemmTypes[indx]));
907db6f9c32SMark Adams #endif
908589dcaf0SStefano Zampini   /* AIR */
909589dcaf0SStefano Zampini #if PETSC_PKG_HYPRE_VERSION_GE(2, 18, 0)
9109566063dSJacob Faibussowitsch   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));
911792fecdfSBarry Smith   PetscCallExternal(HYPRE_BoomerAMGSetRestriction, jac->hsolver, jac->Rtype);
912589dcaf0SStefano Zampini   if (jac->Rtype) {
913589dcaf0SStefano Zampini     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 */
914589dcaf0SStefano Zampini 
9159566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_hypre_boomeramg_strongthresholdR", "Threshold for R", "None", jac->Rstrongthreshold, &jac->Rstrongthreshold, NULL));
916792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetStrongThresholdR, jac->hsolver, jac->Rstrongthreshold);
917589dcaf0SStefano Zampini 
9189566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-pc_hypre_boomeramg_filterthresholdR", "Filter threshold for R", "None", jac->Rfilterthreshold, &jac->Rfilterthreshold, NULL));
919792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetFilterThresholdR, jac->hsolver, jac->Rfilterthreshold);
920589dcaf0SStefano Zampini 
9219566063dSJacob Faibussowitsch     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));
922792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetADropTol, jac->hsolver, jac->Adroptol);
923589dcaf0SStefano Zampini 
9249566063dSJacob Faibussowitsch     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));
925792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetADropType, jac->hsolver, jac->Adroptype);
926589dcaf0SStefano Zampini   }
927589dcaf0SStefano Zampini #endif
928589dcaf0SStefano Zampini 
929ecae95adSPierre Jolivet #if PETSC_PKG_HYPRE_VERSION_LE(9, 9, 9)
93063a3b9bcSJacob Faibussowitsch   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);
931ecae95adSPierre Jolivet #endif
932ecae95adSPierre Jolivet 
9330f1074feSSatish Balay   /* new 3/07 */
934dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_interp_type", "Interpolation type", "None", HYPREBoomerAMGInterpType, PETSC_STATIC_ARRAY_LENGTH(HYPREBoomerAMGInterpType), HYPREBoomerAMGInterpType[0], &indx, &flg));
935589dcaf0SStefano Zampini   if (flg || jac->Rtype) {
936589dcaf0SStefano Zampini     if (flg) jac->interptype = indx;
937792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetInterpType, jac->hsolver, jac->interptype);
9380f1074feSSatish Balay   }
9390f1074feSSatish Balay 
9409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsName("-pc_hypre_boomeramg_print_statistics", "Print statistics", "None", &flg));
94116d9e3a6SLisandro Dalcin   if (flg) {
942b96a4a96SBarry Smith     level = 3;
9439566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_print_statistics", "Print statistics", "None", level, &level, NULL));
9442fa5cd67SKarl Rupp 
945b96a4a96SBarry Smith     jac->printstatistics = PETSC_TRUE;
946792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetPrintLevel, jac->hsolver, level);
9472ae77aedSBarry Smith   }
9482ae77aedSBarry Smith 
9499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsName("-pc_hypre_boomeramg_print_debug", "Print debug information", "None", &flg));
9502ae77aedSBarry Smith   if (flg) {
951b96a4a96SBarry Smith     level = 3;
9529566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_print_debug", "Print debug information", "None", level, &level, NULL));
9532fa5cd67SKarl Rupp 
954b96a4a96SBarry Smith     jac->printstatistics = PETSC_TRUE;
955792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetDebugFlag, jac->hsolver, level);
95616d9e3a6SLisandro Dalcin   }
9578f87f92bSBarry Smith 
9589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg));
9598f87f92bSBarry Smith   if (flg && tmp_truth) {
9608f87f92bSBarry Smith     PetscInt tmp_int;
9619566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", jac->nodal_relax_levels, &tmp_int, &flg));
9628f87f92bSBarry Smith     if (flg) jac->nodal_relax_levels = tmp_int;
963792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetSmoothType, jac->hsolver, 6);
964792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetDomainType, jac->hsolver, 1);
965792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetOverlap, jac->hsolver, 0);
966792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetSmoothNumLevels, jac->hsolver, jac->nodal_relax_levels);
9678f87f92bSBarry Smith   }
9688f87f92bSBarry Smith 
9699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_hypre_boomeramg_keeptranspose", "Avoid transpose matvecs in preconditioner application", "None", jac->keeptranspose, &jac->keeptranspose, NULL));
970792fecdfSBarry Smith   PetscCallExternal(HYPRE_BoomerAMGSetKeepTranspose, jac->hsolver, jac->keeptranspose ? 1 : 0);
971589dcaf0SStefano Zampini 
972589dcaf0SStefano Zampini   /* options for ParaSails solvers */
973dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-pc_hypre_boomeramg_parasails_sym", "Symmetry of matrix and preconditioner", "None", symtlist, PETSC_STATIC_ARRAY_LENGTH(symtlist), symtlist[0], &indx, &flg));
974589dcaf0SStefano Zampini   if (flg) {
975589dcaf0SStefano Zampini     jac->symt = indx;
976792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetSym, jac->hsolver, jac->symt);
977589dcaf0SStefano Zampini   }
978589dcaf0SStefano Zampini 
979d0609cedSBarry Smith   PetscOptionsHeadEnd();
98016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
98116d9e3a6SLisandro Dalcin }
98216d9e3a6SLisandro Dalcin 
9839371c9d4SSatish Balay 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) {
98416d9e3a6SLisandro Dalcin   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
9852cf14000SStefano Zampini   HYPRE_Int oits;
98616d9e3a6SLisandro Dalcin 
98716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9889566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(hypreCitation, &cite));
989792fecdfSBarry Smith   PetscCallExternal(HYPRE_BoomerAMGSetMaxIter, jac->hsolver, its * jac->maxiter);
990792fecdfSBarry Smith   PetscCallExternal(HYPRE_BoomerAMGSetTol, jac->hsolver, rtol);
99116d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_TRUE;
9929566063dSJacob Faibussowitsch   PetscCall(PCApply_HYPRE(pc, b, y));
99316d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_FALSE;
994792fecdfSBarry Smith   PetscCallExternal(HYPRE_BoomerAMGGetNumIterations, jac->hsolver, &oits);
9954d0a8057SBarry Smith   *outits = oits;
9964d0a8057SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
9974d0a8057SBarry Smith   else *reason = PCRICHARDSON_CONVERGED_RTOL;
998792fecdfSBarry Smith   PetscCallExternal(HYPRE_BoomerAMGSetTol, jac->hsolver, jac->tol);
999792fecdfSBarry Smith   PetscCallExternal(HYPRE_BoomerAMGSetMaxIter, jac->hsolver, jac->maxiter);
100016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
100116d9e3a6SLisandro Dalcin }
100216d9e3a6SLisandro Dalcin 
10039371c9d4SSatish Balay static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc, PetscViewer viewer) {
100416d9e3a6SLisandro Dalcin   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
1005ace3abfcSBarry Smith   PetscBool iascii;
100616d9e3a6SLisandro Dalcin 
100716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
10089566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
100916d9e3a6SLisandro Dalcin   if (iascii) {
10109566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  HYPRE BoomerAMG preconditioning\n"));
10119566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Cycle type %s\n", HYPREBoomerAMGCycleType[jac->cycletype]));
101263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Maximum number of levels %" PetscInt_FMT "\n", jac->maxlevels));
101363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Maximum number of iterations PER hypre call %" PetscInt_FMT "\n", jac->maxiter));
10149566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Convergence tolerance PER hypre call %g\n", (double)jac->tol));
10159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Threshold for strong coupling %g\n", (double)jac->strongthreshold));
10169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Interpolation truncation factor %g\n", (double)jac->truncfactor));
101763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Interpolation: max elements per row %" PetscInt_FMT "\n", jac->pmax));
101848a46eb9SPierre Jolivet     if (jac->interp_refine) PetscCall(PetscViewerASCIIPrintf(viewer, "    Interpolation: number of steps of weighted refinement %" PetscInt_FMT "\n", jac->interp_refine));
101963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Number of levels of aggressive coarsening %" PetscInt_FMT "\n", jac->agg_nl));
102063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Number of paths for aggressive coarsening %" PetscInt_FMT "\n", jac->agg_num_paths));
10210f1074feSSatish Balay 
10229566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Maximum row sums %g\n", (double)jac->maxrowsum));
102316d9e3a6SLisandro Dalcin 
102463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Sweeps down         %" PetscInt_FMT "\n", jac->gridsweeps[0]));
102563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Sweeps up           %" PetscInt_FMT "\n", jac->gridsweeps[1]));
102663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Sweeps on coarse    %" PetscInt_FMT "\n", jac->gridsweeps[2]));
102716d9e3a6SLisandro Dalcin 
10289566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Relax down          %s\n", HYPREBoomerAMGRelaxType[jac->relaxtype[0]]));
10299566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Relax up            %s\n", HYPREBoomerAMGRelaxType[jac->relaxtype[1]]));
10309566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Relax on coarse     %s\n", HYPREBoomerAMGRelaxType[jac->relaxtype[2]]));
103116d9e3a6SLisandro Dalcin 
10329566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Relax weight  (all)      %g\n", (double)jac->relaxweight));
10339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Outer relax weight (all) %g\n", (double)jac->outerrelaxweight));
103416d9e3a6SLisandro Dalcin 
103516d9e3a6SLisandro Dalcin     if (jac->relaxorder) {
10369566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Using CF-relaxation\n"));
103716d9e3a6SLisandro Dalcin     } else {
10389566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Not using CF-relaxation\n"));
103916d9e3a6SLisandro Dalcin     }
10406a251517SEike Mueller     if (jac->smoothtype != -1) {
10419566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Smooth type          %s\n", HYPREBoomerAMGSmoothType[jac->smoothtype]));
104263a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Smooth num levels    %" PetscInt_FMT "\n", jac->smoothnumlevels));
10437e352d70SEike Mueller     } else {
10449566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Not using more complex smoothers.\n"));
10451810e44eSEike Mueller     }
10461810e44eSEike Mueller     if (jac->smoothtype == 3) {
104763a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Euclid ILU(k) levels %" PetscInt_FMT "\n", jac->eu_level));
10489566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Euclid ILU(k) drop tolerance %g\n", (double)jac->eu_droptolerance));
104963a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Euclid ILU use Block-Jacobi? %" PetscInt_FMT "\n", jac->eu_bj));
10506a251517SEike Mueller     }
10519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Measure type        %s\n", HYPREBoomerAMGMeasureType[jac->measuretype]));
10529566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Coarsen type        %s\n", HYPREBoomerAMGCoarsenType[jac->coarsentype]));
10539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    Interpolation type  %s\n", jac->interptype != 100 ? HYPREBoomerAMGInterpType[jac->interptype] : "1pt"));
105448a46eb9SPierre Jolivet     if (jac->nodal_coarsening) PetscCall(PetscViewerASCIIPrintf(viewer, "    Using nodal coarsening with HYPRE_BOOMERAMGSetNodal() %" PetscInt_FMT "\n", jac->nodal_coarsening));
10555272c319SBarry Smith     if (jac->vec_interp_variant) {
105663a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    HYPRE_BoomerAMGSetInterpVecVariant() %" PetscInt_FMT "\n", jac->vec_interp_variant));
105763a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    HYPRE_BoomerAMGSetInterpVecQMax() %" PetscInt_FMT "\n", jac->vec_interp_qmax));
10589566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    HYPRE_BoomerAMGSetSmoothInterpVectors() %d\n", jac->vec_interp_smooth));
10598f87f92bSBarry Smith     }
106048a46eb9SPierre Jolivet     if (jac->nodal_relax) PetscCall(PetscViewerASCIIPrintf(viewer, "    Using nodal relaxation via Schwarz smoothing on levels %" PetscInt_FMT "\n", jac->nodal_relax_levels));
1061db6f9c32SMark Adams #if PETSC_PKG_HYPRE_VERSION_GE(2, 23, 0)
10629566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    SpGEMM type         %s\n", jac->spgemm_type));
1063db6f9c32SMark Adams #endif
1064589dcaf0SStefano Zampini     /* AIR */
1065589dcaf0SStefano Zampini     if (jac->Rtype) {
106663a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Using approximate ideal restriction type %" PetscInt_FMT "\n", jac->Rtype));
10679566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "      Threshold for R %g\n", (double)jac->Rstrongthreshold));
10689566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "      Filter for R %g\n", (double)jac->Rfilterthreshold));
10699566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "      A drop tolerance %g\n", (double)jac->Adroptol));
107063a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "      A drop type %" PetscInt_FMT "\n", jac->Adroptype));
1071589dcaf0SStefano Zampini     }
107216d9e3a6SLisandro Dalcin   }
107316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
107416d9e3a6SLisandro Dalcin }
107516d9e3a6SLisandro Dalcin 
10769371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc, PetscOptionItems *PetscOptionsObject) {
107716d9e3a6SLisandro Dalcin   PC_HYPRE   *jac = (PC_HYPRE *)pc->data;
10784ddd07fcSJed Brown   PetscInt    indx;
1079ace3abfcSBarry Smith   PetscBool   flag;
108016d9e3a6SLisandro Dalcin   const char *symtlist[] = {"nonsymmetric", "SPD", "nonsymmetric,SPD"};
108116d9e3a6SLisandro Dalcin 
108216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
1083d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "HYPRE ParaSails Options");
10849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_parasails_nlevels", "Number of number of levels", "None", jac->nlevels, &jac->nlevels, 0));
10859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_parasails_thresh", "Threshold", "None", jac->threshold, &jac->threshold, &flag));
1086792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_ParaSailsSetParams, jac->hsolver, jac->threshold, jac->nlevels);
108716d9e3a6SLisandro Dalcin 
10889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_parasails_filter", "filter", "None", jac->filter, &jac->filter, &flag));
1089792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_ParaSailsSetFilter, jac->hsolver, jac->filter);
109016d9e3a6SLisandro Dalcin 
10919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_parasails_loadbal", "Load balance", "None", jac->loadbal, &jac->loadbal, &flag));
1092792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_ParaSailsSetLoadbal, jac->hsolver, jac->loadbal);
109316d9e3a6SLisandro Dalcin 
10949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_hypre_parasails_logging", "Print info to screen", "None", (PetscBool)jac->logging, (PetscBool *)&jac->logging, &flag));
1095792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_ParaSailsSetLogging, jac->hsolver, jac->logging);
109616d9e3a6SLisandro Dalcin 
10979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_hypre_parasails_reuse", "Reuse nonzero pattern in preconditioner", "None", (PetscBool)jac->ruse, (PetscBool *)&jac->ruse, &flag));
1098792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_ParaSailsSetReuse, jac->hsolver, jac->ruse);
109916d9e3a6SLisandro Dalcin 
1100dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-pc_hypre_parasails_sym", "Symmetry of matrix and preconditioner", "None", symtlist, PETSC_STATIC_ARRAY_LENGTH(symtlist), symtlist[0], &indx, &flag));
110116d9e3a6SLisandro Dalcin   if (flag) {
110216d9e3a6SLisandro Dalcin     jac->symt = indx;
1103792fecdfSBarry Smith     PetscCallExternal(HYPRE_ParaSailsSetSym, jac->hsolver, jac->symt);
110416d9e3a6SLisandro Dalcin   }
110516d9e3a6SLisandro Dalcin 
1106d0609cedSBarry Smith   PetscOptionsHeadEnd();
110716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
110816d9e3a6SLisandro Dalcin }
110916d9e3a6SLisandro Dalcin 
11109371c9d4SSatish Balay static PetscErrorCode PCView_HYPRE_ParaSails(PC pc, PetscViewer viewer) {
111116d9e3a6SLisandro Dalcin   PC_HYPRE   *jac = (PC_HYPRE *)pc->data;
1112ace3abfcSBarry Smith   PetscBool   iascii;
1113feb237baSPierre Jolivet   const char *symt = 0;
111416d9e3a6SLisandro Dalcin 
111516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
11169566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
111716d9e3a6SLisandro Dalcin   if (iascii) {
11189566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  HYPRE ParaSails preconditioning\n"));
111963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    nlevels %" PetscInt_FMT "\n", jac->nlevels));
11209566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    threshold %g\n", (double)jac->threshold));
11219566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    filter %g\n", (double)jac->filter));
11229566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    load balance %g\n", (double)jac->loadbal));
11239566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    reuse nonzero structure %s\n", PetscBools[jac->ruse]));
11249566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    print info to screen %s\n", PetscBools[jac->logging]));
11252fa5cd67SKarl Rupp     if (!jac->symt) symt = "nonsymmetric matrix and preconditioner";
11262fa5cd67SKarl Rupp     else if (jac->symt == 1) symt = "SPD matrix and preconditioner";
11272fa5cd67SKarl Rupp     else if (jac->symt == 2) symt = "nonsymmetric matrix but SPD preconditioner";
112863a3b9bcSJacob Faibussowitsch     else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Unknown HYPRE ParaSails symmetric option %" PetscInt_FMT, jac->symt);
11299566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    %s\n", symt));
113016d9e3a6SLisandro Dalcin   }
113116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
113216d9e3a6SLisandro Dalcin }
1133f1580f4eSBarry Smith 
11349371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_HYPRE_AMS(PC pc, PetscOptionItems *PetscOptionsObject) {
11354cb006feSStefano Zampini   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
11364cb006feSStefano Zampini   PetscInt  n;
11374cb006feSStefano Zampini   PetscBool flag, flag2, flag3, flag4;
11384cb006feSStefano Zampini 
11394cb006feSStefano Zampini   PetscFunctionBegin;
1140d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "HYPRE AMS Options");
11419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_ams_print_level", "Debugging output level for AMS", "None", jac->as_print, &jac->as_print, &flag));
1142792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_AMSSetPrintLevel, jac->hsolver, jac->as_print);
11439566063dSJacob Faibussowitsch   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));
1144792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_AMSSetMaxIter, jac->hsolver, jac->as_max_iter);
11459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_ams_cycle_type", "Cycle type for AMS multigrid", "None", jac->ams_cycle_type, &jac->ams_cycle_type, &flag));
1146792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_AMSSetCycleType, jac->hsolver, jac->ams_cycle_type);
11479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_ams_tol", "Error tolerance for AMS multigrid", "None", jac->as_tol, &jac->as_tol, &flag));
1148792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_AMSSetTol, jac->hsolver, jac->as_tol);
11499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_ams_relax_type", "Relaxation type for AMS smoother", "None", jac->as_relax_type, &jac->as_relax_type, &flag));
11509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_ams_relax_times", "Number of relaxation steps for AMS smoother", "None", jac->as_relax_times, &jac->as_relax_times, &flag2));
11519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_ams_relax_weight", "Relaxation weight for AMS smoother", "None", jac->as_relax_weight, &jac->as_relax_weight, &flag3));
11529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_ams_omega", "SSOR coefficient for AMS smoother", "None", jac->as_omega, &jac->as_omega, &flag4));
115348a46eb9SPierre Jolivet   if (flag || flag2 || flag3 || flag4) PetscCallExternal(HYPRE_AMSSetSmoothingOptions, jac->hsolver, jac->as_relax_type, jac->as_relax_times, jac->as_relax_weight, jac->as_omega);
11549566063dSJacob Faibussowitsch   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));
11554cb006feSStefano Zampini   n = 5;
11569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-pc_hypre_ams_amg_alpha_options", "AMG options for vector Poisson", "None", jac->as_amg_alpha_opts, &n, &flag2));
11574cb006feSStefano Zampini   if (flag || flag2) {
1158792fecdfSBarry Smith     PetscCallExternal(HYPRE_AMSSetAlphaAMGOptions, jac->hsolver, jac->as_amg_alpha_opts[0], /* AMG coarsen type */
1159863406b8SStefano Zampini                       jac->as_amg_alpha_opts[1],                                            /* AMG agg_levels */
1160863406b8SStefano Zampini                       jac->as_amg_alpha_opts[2],                                            /* AMG relax_type */
11619371c9d4SSatish Balay                       jac->as_amg_alpha_theta, jac->as_amg_alpha_opts[3],                   /* AMG interp_type */
1162a74df02fSJacob Faibussowitsch                       jac->as_amg_alpha_opts[4]);                                           /* AMG Pmax */
11634cb006feSStefano Zampini   }
11649566063dSJacob Faibussowitsch   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));
11654cb006feSStefano Zampini   n = 5;
11669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-pc_hypre_ams_amg_beta_options", "AMG options for scalar Poisson solver", "None", jac->as_amg_beta_opts, &n, &flag2));
11674cb006feSStefano Zampini   if (flag || flag2) {
1168792fecdfSBarry Smith     PetscCallExternal(HYPRE_AMSSetBetaAMGOptions, jac->hsolver, jac->as_amg_beta_opts[0], /* AMG coarsen type */
1169863406b8SStefano Zampini                       jac->as_amg_beta_opts[1],                                           /* AMG agg_levels */
1170863406b8SStefano Zampini                       jac->as_amg_beta_opts[2],                                           /* AMG relax_type */
11719371c9d4SSatish Balay                       jac->as_amg_beta_theta, jac->as_amg_beta_opts[3],                   /* AMG interp_type */
1172a74df02fSJacob Faibussowitsch                       jac->as_amg_beta_opts[4]);                                          /* AMG Pmax */
11734cb006feSStefano Zampini   }
11749566063dSJacob Faibussowitsch   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));
117523df4f25SStefano Zampini   if (flag) { /* override HYPRE's default only if the options is used */
1176792fecdfSBarry Smith     PetscCallExternal(HYPRE_AMSSetProjectionFrequency, jac->hsolver, jac->ams_proj_freq);
117723df4f25SStefano Zampini   }
1178d0609cedSBarry Smith   PetscOptionsHeadEnd();
11794cb006feSStefano Zampini   PetscFunctionReturn(0);
11804cb006feSStefano Zampini }
11814cb006feSStefano Zampini 
11829371c9d4SSatish Balay static PetscErrorCode PCView_HYPRE_AMS(PC pc, PetscViewer viewer) {
11834cb006feSStefano Zampini   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
11844cb006feSStefano Zampini   PetscBool iascii;
11854cb006feSStefano Zampini 
11864cb006feSStefano Zampini   PetscFunctionBegin;
11879566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
11884cb006feSStefano Zampini   if (iascii) {
11899566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  HYPRE AMS preconditioning\n"));
119063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    subspace iterations per application %" PetscInt_FMT "\n", jac->as_max_iter));
119163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    subspace cycle type %" PetscInt_FMT "\n", jac->ams_cycle_type));
119263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    subspace iteration tolerance %g\n", (double)jac->as_tol));
119363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    smoother type %" PetscInt_FMT "\n", jac->as_relax_type));
119463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    number of smoothing steps %" PetscInt_FMT "\n", jac->as_relax_times));
119563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    smoother weight %g\n", (double)jac->as_relax_weight));
119663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    smoother omega %g\n", (double)jac->as_omega));
11974cb006feSStefano Zampini     if (jac->alpha_Poisson) {
11989566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    vector Poisson solver (passed in by user)\n"));
11994cb006feSStefano Zampini     } else {
12009566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    vector Poisson solver (computed) \n"));
12014cb006feSStefano Zampini     }
120263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        boomerAMG coarsening type %" PetscInt_FMT "\n", jac->as_amg_alpha_opts[0]));
120363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        boomerAMG levels of aggressive coarsening %" PetscInt_FMT "\n", jac->as_amg_alpha_opts[1]));
120463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        boomerAMG relaxation type %" PetscInt_FMT "\n", jac->as_amg_alpha_opts[2]));
120563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        boomerAMG interpolation type %" PetscInt_FMT "\n", jac->as_amg_alpha_opts[3]));
120663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        boomerAMG max nonzero elements in interpolation rows %" PetscInt_FMT "\n", jac->as_amg_alpha_opts[4]));
120763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        boomerAMG strength threshold %g\n", (double)jac->as_amg_alpha_theta));
12084cb006feSStefano Zampini     if (!jac->ams_beta_is_zero) {
12094cb006feSStefano Zampini       if (jac->beta_Poisson) {
12109566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "    scalar Poisson solver (passed in by user)\n"));
12114cb006feSStefano Zampini       } else {
12129566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "    scalar Poisson solver (computed) \n"));
12134cb006feSStefano Zampini       }
121463a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "        boomerAMG coarsening type %" PetscInt_FMT "\n", jac->as_amg_beta_opts[0]));
121563a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "        boomerAMG levels of aggressive coarsening %" PetscInt_FMT "\n", jac->as_amg_beta_opts[1]));
121663a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "        boomerAMG relaxation type %" PetscInt_FMT "\n", jac->as_amg_beta_opts[2]));
121763a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "        boomerAMG interpolation type %" PetscInt_FMT "\n", jac->as_amg_beta_opts[3]));
121863a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "        boomerAMG max nonzero elements in interpolation rows %" PetscInt_FMT "\n", jac->as_amg_beta_opts[4]));
121963a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "        boomerAMG strength threshold %g\n", (double)jac->as_amg_beta_theta));
122048a46eb9SPierre Jolivet       if (jac->ams_beta_is_zero_part) PetscCall(PetscViewerASCIIPrintf(viewer, "        compatible subspace projection frequency %" PetscInt_FMT " (-1 HYPRE uses default)\n", jac->ams_proj_freq));
122123df4f25SStefano Zampini     } else {
12229566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    scalar Poisson solver not used (zero-conductivity everywhere) \n"));
12234cb006feSStefano Zampini     }
12244cb006feSStefano Zampini   }
12254cb006feSStefano Zampini   PetscFunctionReturn(0);
12264cb006feSStefano Zampini }
12274cb006feSStefano Zampini 
12289371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_HYPRE_ADS(PC pc, PetscOptionItems *PetscOptionsObject) {
1229863406b8SStefano Zampini   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
1230863406b8SStefano Zampini   PetscInt  n;
1231863406b8SStefano Zampini   PetscBool flag, flag2, flag3, flag4;
1232863406b8SStefano Zampini 
1233863406b8SStefano Zampini   PetscFunctionBegin;
1234d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "HYPRE ADS Options");
12359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_ads_print_level", "Debugging output level for ADS", "None", jac->as_print, &jac->as_print, &flag));
1236792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_ADSSetPrintLevel, jac->hsolver, jac->as_print);
12379566063dSJacob Faibussowitsch   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));
1238792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_ADSSetMaxIter, jac->hsolver, jac->as_max_iter);
12399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_ads_cycle_type", "Cycle type for ADS multigrid", "None", jac->ads_cycle_type, &jac->ads_cycle_type, &flag));
1240792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_ADSSetCycleType, jac->hsolver, jac->ads_cycle_type);
12419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_ads_tol", "Error tolerance for ADS multigrid", "None", jac->as_tol, &jac->as_tol, &flag));
1242792fecdfSBarry Smith   if (flag) PetscCallExternal(HYPRE_ADSSetTol, jac->hsolver, jac->as_tol);
12439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_ads_relax_type", "Relaxation type for ADS smoother", "None", jac->as_relax_type, &jac->as_relax_type, &flag));
12449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hypre_ads_relax_times", "Number of relaxation steps for ADS smoother", "None", jac->as_relax_times, &jac->as_relax_times, &flag2));
12459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_ads_relax_weight", "Relaxation weight for ADS smoother", "None", jac->as_relax_weight, &jac->as_relax_weight, &flag3));
12469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_hypre_ads_omega", "SSOR coefficient for ADS smoother", "None", jac->as_omega, &jac->as_omega, &flag4));
124748a46eb9SPierre Jolivet   if (flag || flag2 || flag3 || flag4) PetscCallExternal(HYPRE_ADSSetSmoothingOptions, jac->hsolver, jac->as_relax_type, jac->as_relax_times, jac->as_relax_weight, jac->as_omega);
12489566063dSJacob Faibussowitsch   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));
1249863406b8SStefano Zampini   n = 5;
12509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-pc_hypre_ads_ams_options", "AMG options for AMS solver inside ADS", "None", jac->as_amg_alpha_opts, &n, &flag2));
12519566063dSJacob Faibussowitsch   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));
1252863406b8SStefano Zampini   if (flag || flag2 || flag3) {
1253792fecdfSBarry Smith     PetscCallExternal(HYPRE_ADSSetAMSOptions, jac->hsolver, jac->ams_cycle_type, /* AMS cycle type */
1254863406b8SStefano Zampini                       jac->as_amg_alpha_opts[0],                                 /* AMG coarsen type */
1255863406b8SStefano Zampini                       jac->as_amg_alpha_opts[1],                                 /* AMG agg_levels */
1256863406b8SStefano Zampini                       jac->as_amg_alpha_opts[2],                                 /* AMG relax_type */
12579371c9d4SSatish Balay                       jac->as_amg_alpha_theta, jac->as_amg_alpha_opts[3],        /* AMG interp_type */
1258a74df02fSJacob Faibussowitsch                       jac->as_amg_alpha_opts[4]);                                /* AMG Pmax */
1259863406b8SStefano Zampini   }
12609566063dSJacob Faibussowitsch   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));
1261863406b8SStefano Zampini   n = 5;
12629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-pc_hypre_ads_amg_options", "AMG options for vector AMG solver inside ADS", "None", jac->as_amg_beta_opts, &n, &flag2));
1263863406b8SStefano Zampini   if (flag || flag2) {
1264792fecdfSBarry Smith     PetscCallExternal(HYPRE_ADSSetAMGOptions, jac->hsolver, jac->as_amg_beta_opts[0], /* AMG coarsen type */
1265863406b8SStefano Zampini                       jac->as_amg_beta_opts[1],                                       /* AMG agg_levels */
1266863406b8SStefano Zampini                       jac->as_amg_beta_opts[2],                                       /* AMG relax_type */
12679371c9d4SSatish Balay                       jac->as_amg_beta_theta, jac->as_amg_beta_opts[3],               /* AMG interp_type */
1268a74df02fSJacob Faibussowitsch                       jac->as_amg_beta_opts[4]);                                      /* AMG Pmax */
1269863406b8SStefano Zampini   }
1270d0609cedSBarry Smith   PetscOptionsHeadEnd();
1271863406b8SStefano Zampini   PetscFunctionReturn(0);
1272863406b8SStefano Zampini }
1273863406b8SStefano Zampini 
12749371c9d4SSatish Balay static PetscErrorCode PCView_HYPRE_ADS(PC pc, PetscViewer viewer) {
1275863406b8SStefano Zampini   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
1276863406b8SStefano Zampini   PetscBool iascii;
1277863406b8SStefano Zampini 
1278863406b8SStefano Zampini   PetscFunctionBegin;
12799566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1280863406b8SStefano Zampini   if (iascii) {
12819566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  HYPRE ADS preconditioning\n"));
128263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    subspace iterations per application %" PetscInt_FMT "\n", jac->as_max_iter));
128363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    subspace cycle type %" PetscInt_FMT "\n", jac->ads_cycle_type));
128463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    subspace iteration tolerance %g\n", (double)jac->as_tol));
128563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    smoother type %" PetscInt_FMT "\n", jac->as_relax_type));
128663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    number of smoothing steps %" PetscInt_FMT "\n", jac->as_relax_times));
128763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    smoother weight %g\n", (double)jac->as_relax_weight));
128863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    smoother omega %g\n", (double)jac->as_omega));
12899566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    AMS solver using boomerAMG\n"));
129063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        subspace cycle type %" PetscInt_FMT "\n", jac->ams_cycle_type));
129163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        coarsening type %" PetscInt_FMT "\n", jac->as_amg_alpha_opts[0]));
129263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        levels of aggressive coarsening %" PetscInt_FMT "\n", jac->as_amg_alpha_opts[1]));
129363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        relaxation type %" PetscInt_FMT "\n", jac->as_amg_alpha_opts[2]));
129463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        interpolation type %" PetscInt_FMT "\n", jac->as_amg_alpha_opts[3]));
129563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        max nonzero elements in interpolation rows %" PetscInt_FMT "\n", jac->as_amg_alpha_opts[4]));
129663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        strength threshold %g\n", (double)jac->as_amg_alpha_theta));
12979566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    vector Poisson solver using boomerAMG\n"));
129863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        coarsening type %" PetscInt_FMT "\n", jac->as_amg_beta_opts[0]));
129963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        levels of aggressive coarsening %" PetscInt_FMT "\n", jac->as_amg_beta_opts[1]));
130063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        relaxation type %" PetscInt_FMT "\n", jac->as_amg_beta_opts[2]));
130163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        interpolation type %" PetscInt_FMT "\n", jac->as_amg_beta_opts[3]));
130263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        max nonzero elements in interpolation rows %" PetscInt_FMT "\n", jac->as_amg_beta_opts[4]));
130363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "        strength threshold %g\n", (double)jac->as_amg_beta_theta));
1304863406b8SStefano Zampini   }
1305863406b8SStefano Zampini   PetscFunctionReturn(0);
1306863406b8SStefano Zampini }
1307863406b8SStefano Zampini 
13089371c9d4SSatish Balay static PetscErrorCode PCHYPRESetDiscreteGradient_HYPRE(PC pc, Mat G) {
13094cb006feSStefano Zampini   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
13105ac14e1cSStefano Zampini   PetscBool ishypre;
13114cb006feSStefano Zampini 
13124cb006feSStefano Zampini   PetscFunctionBegin;
13139566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)G, MATHYPRE, &ishypre));
13145ac14e1cSStefano Zampini   if (ishypre) {
13159566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)G));
13169566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&jac->G));
13175ac14e1cSStefano Zampini     jac->G = G;
13185ac14e1cSStefano Zampini   } else {
13199566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&jac->G));
13209566063dSJacob Faibussowitsch     PetscCall(MatConvert(G, MATHYPRE, MAT_INITIAL_MATRIX, &jac->G));
13215ac14e1cSStefano Zampini   }
13224cb006feSStefano Zampini   PetscFunctionReturn(0);
13234cb006feSStefano Zampini }
13244cb006feSStefano Zampini 
13254cb006feSStefano Zampini /*@
1326f1580f4eSBarry Smith    PCHYPRESetDiscreteGradient - Set discrete gradient matrix for `PCHYPRE` type of ams or ads
13274cb006feSStefano Zampini 
1328f1580f4eSBarry Smith    Collective on pc
13294cb006feSStefano Zampini 
13304cb006feSStefano Zampini    Input Parameters:
13314cb006feSStefano Zampini +  pc - the preconditioning context
13324cb006feSStefano Zampini -  G - the discrete gradient
13334cb006feSStefano Zampini 
13344cb006feSStefano Zampini    Level: intermediate
13354cb006feSStefano Zampini 
133695452b02SPatrick Sanan    Notes:
133795452b02SPatrick Sanan     G should have as many rows as the number of edges and as many columns as the number of vertices in the mesh
1338147403d9SBarry Smith 
1339863406b8SStefano Zampini     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
13404cb006feSStefano Zampini 
1341f1580f4eSBarry Smith    Developer Note:
1342f1580f4eSBarry Smith    This automatically converts the matrix to `MATHYPRE` if it is not already of that type
1343f1580f4eSBarry Smith 
1344f1580f4eSBarry Smith .seealso: `PCHYPRE`, `PCHYPRESetDiscreteCurl()`
13454cb006feSStefano Zampini @*/
13469371c9d4SSatish Balay PetscErrorCode PCHYPRESetDiscreteGradient(PC pc, Mat G) {
13474cb006feSStefano Zampini   PetscFunctionBegin;
13484cb006feSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
13494cb006feSStefano Zampini   PetscValidHeaderSpecific(G, MAT_CLASSID, 2);
13504cb006feSStefano Zampini   PetscCheckSameComm(pc, 1, G, 2);
1351cac4c232SBarry Smith   PetscTryMethod(pc, "PCHYPRESetDiscreteGradient_C", (PC, Mat), (pc, G));
13524cb006feSStefano Zampini   PetscFunctionReturn(0);
13534cb006feSStefano Zampini }
13544cb006feSStefano Zampini 
13559371c9d4SSatish Balay static PetscErrorCode PCHYPRESetDiscreteCurl_HYPRE(PC pc, Mat C) {
1356863406b8SStefano Zampini   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
13575ac14e1cSStefano Zampini   PetscBool ishypre;
1358863406b8SStefano Zampini 
1359863406b8SStefano Zampini   PetscFunctionBegin;
13609566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)C, MATHYPRE, &ishypre));
13615ac14e1cSStefano Zampini   if (ishypre) {
13629566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)C));
13639566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&jac->C));
13645ac14e1cSStefano Zampini     jac->C = C;
13655ac14e1cSStefano Zampini   } else {
13669566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&jac->C));
13679566063dSJacob Faibussowitsch     PetscCall(MatConvert(C, MATHYPRE, MAT_INITIAL_MATRIX, &jac->C));
13685ac14e1cSStefano Zampini   }
1369863406b8SStefano Zampini   PetscFunctionReturn(0);
1370863406b8SStefano Zampini }
1371863406b8SStefano Zampini 
1372863406b8SStefano Zampini /*@
1373f1580f4eSBarry Smith    PCHYPRESetDiscreteCurl - Set discrete curl matrx for `PCHYPRE` type of ads
1374863406b8SStefano Zampini 
1375f1580f4eSBarry Smith    Collective on pc
1376863406b8SStefano Zampini 
1377863406b8SStefano Zampini    Input Parameters:
1378863406b8SStefano Zampini +  pc - the preconditioning context
1379863406b8SStefano Zampini -  C - the discrete curl
1380863406b8SStefano Zampini 
1381863406b8SStefano Zampini    Level: intermediate
1382863406b8SStefano Zampini 
138395452b02SPatrick Sanan    Notes:
138495452b02SPatrick Sanan     C should have as many rows as the number of faces and as many columns as the number of edges in the mesh
1385147403d9SBarry Smith 
1386863406b8SStefano Zampini     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
1387863406b8SStefano Zampini 
1388f1580f4eSBarry Smith    Developer Note:
1389f1580f4eSBarry Smith    This automatically converts the matrix to `MATHYPRE` if it is not already of that type
1390f1580f4eSBarry Smith 
1391f1580f4eSBarry Smith    If this is only for  `PCHYPRE` type of ads it should be called `PCHYPREADSSetDiscreteCurl()`
1392f1580f4eSBarry Smith 
1393f1580f4eSBarry Smith .seealso: `PCHYPRE`, `PCHYPRESetDiscreteGradient()`
1394863406b8SStefano Zampini @*/
13959371c9d4SSatish Balay PetscErrorCode PCHYPRESetDiscreteCurl(PC pc, Mat C) {
1396863406b8SStefano Zampini   PetscFunctionBegin;
1397863406b8SStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1398863406b8SStefano Zampini   PetscValidHeaderSpecific(C, MAT_CLASSID, 2);
1399863406b8SStefano Zampini   PetscCheckSameComm(pc, 1, C, 2);
1400cac4c232SBarry Smith   PetscTryMethod(pc, "PCHYPRESetDiscreteCurl_C", (PC, Mat), (pc, C));
1401863406b8SStefano Zampini   PetscFunctionReturn(0);
1402863406b8SStefano Zampini }
1403863406b8SStefano Zampini 
14049371c9d4SSatish Balay static PetscErrorCode PCHYPRESetInterpolations_HYPRE(PC pc, PetscInt dim, Mat RT_PiFull, Mat RT_Pi[], Mat ND_PiFull, Mat ND_Pi[]) {
14056bf688a0SCe Qin   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
14066bf688a0SCe Qin   PetscBool ishypre;
14076bf688a0SCe Qin   PetscInt  i;
14086bf688a0SCe Qin   PetscFunctionBegin;
14096bf688a0SCe Qin 
14109566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->RT_PiFull));
14119566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->ND_PiFull));
14126bf688a0SCe Qin   for (i = 0; i < 3; ++i) {
14139566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&jac->RT_Pi[i]));
14149566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&jac->ND_Pi[i]));
14156bf688a0SCe Qin   }
14166bf688a0SCe Qin 
14176bf688a0SCe Qin   jac->dim = dim;
14186bf688a0SCe Qin   if (RT_PiFull) {
14199566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)RT_PiFull, MATHYPRE, &ishypre));
14206bf688a0SCe Qin     if (ishypre) {
14219566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)RT_PiFull));
14226bf688a0SCe Qin       jac->RT_PiFull = RT_PiFull;
14236bf688a0SCe Qin     } else {
14249566063dSJacob Faibussowitsch       PetscCall(MatConvert(RT_PiFull, MATHYPRE, MAT_INITIAL_MATRIX, &jac->RT_PiFull));
14256bf688a0SCe Qin     }
14266bf688a0SCe Qin   }
14276bf688a0SCe Qin   if (RT_Pi) {
14286bf688a0SCe Qin     for (i = 0; i < dim; ++i) {
14296bf688a0SCe Qin       if (RT_Pi[i]) {
14309566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)RT_Pi[i], MATHYPRE, &ishypre));
14316bf688a0SCe Qin         if (ishypre) {
14329566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)RT_Pi[i]));
14336bf688a0SCe Qin           jac->RT_Pi[i] = RT_Pi[i];
14346bf688a0SCe Qin         } else {
14359566063dSJacob Faibussowitsch           PetscCall(MatConvert(RT_Pi[i], MATHYPRE, MAT_INITIAL_MATRIX, &jac->RT_Pi[i]));
14366bf688a0SCe Qin         }
14376bf688a0SCe Qin       }
14386bf688a0SCe Qin     }
14396bf688a0SCe Qin   }
14406bf688a0SCe Qin   if (ND_PiFull) {
14419566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)ND_PiFull, MATHYPRE, &ishypre));
14426bf688a0SCe Qin     if (ishypre) {
14439566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)ND_PiFull));
14446bf688a0SCe Qin       jac->ND_PiFull = ND_PiFull;
14456bf688a0SCe Qin     } else {
14469566063dSJacob Faibussowitsch       PetscCall(MatConvert(ND_PiFull, MATHYPRE, MAT_INITIAL_MATRIX, &jac->ND_PiFull));
14476bf688a0SCe Qin     }
14486bf688a0SCe Qin   }
14496bf688a0SCe Qin   if (ND_Pi) {
14506bf688a0SCe Qin     for (i = 0; i < dim; ++i) {
14516bf688a0SCe Qin       if (ND_Pi[i]) {
14529566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)ND_Pi[i], MATHYPRE, &ishypre));
14536bf688a0SCe Qin         if (ishypre) {
14549566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)ND_Pi[i]));
14556bf688a0SCe Qin           jac->ND_Pi[i] = ND_Pi[i];
14566bf688a0SCe Qin         } else {
14579566063dSJacob Faibussowitsch           PetscCall(MatConvert(ND_Pi[i], MATHYPRE, MAT_INITIAL_MATRIX, &jac->ND_Pi[i]));
14586bf688a0SCe Qin         }
14596bf688a0SCe Qin       }
14606bf688a0SCe Qin     }
14616bf688a0SCe Qin   }
14626bf688a0SCe Qin 
14636bf688a0SCe Qin   PetscFunctionReturn(0);
14646bf688a0SCe Qin }
14656bf688a0SCe Qin 
14666bf688a0SCe Qin /*@
1467f1580f4eSBarry Smith    PCHYPRESetInterpolations - Set interpolation matrices for `PCHYPRE` type of ams or ads
14686bf688a0SCe Qin 
1469f1580f4eSBarry Smith    Collective on pc
14706bf688a0SCe Qin 
14716bf688a0SCe Qin    Input Parameters:
14726bf688a0SCe Qin +  pc - the preconditioning context
14736bf688a0SCe Qin -  dim - the dimension of the problem, only used in AMS
14746bf688a0SCe Qin -  RT_PiFull - Raviart-Thomas interpolation matrix
14756bf688a0SCe Qin -  RT_Pi - x/y/z component of Raviart-Thomas interpolation matrix
14766bf688a0SCe Qin -  ND_PiFull - Nedelec interpolation matrix
14776bf688a0SCe Qin -  ND_Pi - x/y/z component of Nedelec interpolation matrix
14786bf688a0SCe Qin 
1479f1580f4eSBarry Smith    Level: intermediate
1480f1580f4eSBarry Smith 
148195452b02SPatrick Sanan    Notes:
148295452b02SPatrick Sanan     For AMS, only Nedelec interpolation matrices are needed, the Raviart-Thomas interpolation matrices can be set to NULL.
1483147403d9SBarry Smith 
14846bf688a0SCe Qin     For ADS, both type of interpolation matrices are needed.
1485147403d9SBarry Smith 
1486f1580f4eSBarry Smith    Developer Note:
1487f1580f4eSBarry Smith    This automatically converts the matrix to `MATHYPRE` if it is not already of that type
14886bf688a0SCe Qin 
1489f1580f4eSBarry Smith .seealso: `PCHYPRE`
14906bf688a0SCe Qin @*/
14919371c9d4SSatish Balay PetscErrorCode PCHYPRESetInterpolations(PC pc, PetscInt dim, Mat RT_PiFull, Mat RT_Pi[], Mat ND_PiFull, Mat ND_Pi[]) {
14926bf688a0SCe Qin   PetscInt i;
14936bf688a0SCe Qin 
14946bf688a0SCe Qin   PetscFunctionBegin;
14956bf688a0SCe Qin   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
14966bf688a0SCe Qin   if (RT_PiFull) {
14976bf688a0SCe Qin     PetscValidHeaderSpecific(RT_PiFull, MAT_CLASSID, 3);
14986bf688a0SCe Qin     PetscCheckSameComm(pc, 1, RT_PiFull, 3);
14996bf688a0SCe Qin   }
15006bf688a0SCe Qin   if (RT_Pi) {
15016bf688a0SCe Qin     PetscValidPointer(RT_Pi, 4);
15026bf688a0SCe Qin     for (i = 0; i < dim; ++i) {
15036bf688a0SCe Qin       if (RT_Pi[i]) {
15046bf688a0SCe Qin         PetscValidHeaderSpecific(RT_Pi[i], MAT_CLASSID, 4);
15056bf688a0SCe Qin         PetscCheckSameComm(pc, 1, RT_Pi[i], 4);
15066bf688a0SCe Qin       }
15076bf688a0SCe Qin     }
15086bf688a0SCe Qin   }
15096bf688a0SCe Qin   if (ND_PiFull) {
15106bf688a0SCe Qin     PetscValidHeaderSpecific(ND_PiFull, MAT_CLASSID, 5);
15116bf688a0SCe Qin     PetscCheckSameComm(pc, 1, ND_PiFull, 5);
15126bf688a0SCe Qin   }
15136bf688a0SCe Qin   if (ND_Pi) {
15146bf688a0SCe Qin     PetscValidPointer(ND_Pi, 6);
15156bf688a0SCe Qin     for (i = 0; i < dim; ++i) {
15166bf688a0SCe Qin       if (ND_Pi[i]) {
15176bf688a0SCe Qin         PetscValidHeaderSpecific(ND_Pi[i], MAT_CLASSID, 6);
15186bf688a0SCe Qin         PetscCheckSameComm(pc, 1, ND_Pi[i], 6);
15196bf688a0SCe Qin       }
15206bf688a0SCe Qin     }
15216bf688a0SCe Qin   }
1522cac4c232SBarry Smith   PetscTryMethod(pc, "PCHYPRESetInterpolations_C", (PC, PetscInt, Mat, Mat[], Mat, Mat[]), (pc, dim, RT_PiFull, RT_Pi, ND_PiFull, ND_Pi));
15236bf688a0SCe Qin   PetscFunctionReturn(0);
15246bf688a0SCe Qin }
15256bf688a0SCe Qin 
15269371c9d4SSatish Balay static PetscErrorCode PCHYPRESetPoissonMatrix_HYPRE(PC pc, Mat A, PetscBool isalpha) {
15274cb006feSStefano Zampini   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
15285ac14e1cSStefano Zampini   PetscBool ishypre;
15294cb006feSStefano Zampini 
15304cb006feSStefano Zampini   PetscFunctionBegin;
15319566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATHYPRE, &ishypre));
15325ac14e1cSStefano Zampini   if (ishypre) {
15335ac14e1cSStefano Zampini     if (isalpha) {
15349566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)A));
15359566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&jac->alpha_Poisson));
15365ac14e1cSStefano Zampini       jac->alpha_Poisson = A;
15375ac14e1cSStefano Zampini     } else {
15385ac14e1cSStefano Zampini       if (A) {
15399566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)A));
15405ac14e1cSStefano Zampini       } else {
15415ac14e1cSStefano Zampini         jac->ams_beta_is_zero = PETSC_TRUE;
15425ac14e1cSStefano Zampini       }
15439566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&jac->beta_Poisson));
15445ac14e1cSStefano Zampini       jac->beta_Poisson = A;
15455ac14e1cSStefano Zampini     }
15465ac14e1cSStefano Zampini   } else {
15475ac14e1cSStefano Zampini     if (isalpha) {
15489566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&jac->alpha_Poisson));
15499566063dSJacob Faibussowitsch       PetscCall(MatConvert(A, MATHYPRE, MAT_INITIAL_MATRIX, &jac->alpha_Poisson));
15505ac14e1cSStefano Zampini     } else {
15515ac14e1cSStefano Zampini       if (A) {
15529566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&jac->beta_Poisson));
15539566063dSJacob Faibussowitsch         PetscCall(MatConvert(A, MATHYPRE, MAT_INITIAL_MATRIX, &jac->beta_Poisson));
15545ac14e1cSStefano Zampini       } else {
15559566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&jac->beta_Poisson));
15565ac14e1cSStefano Zampini         jac->ams_beta_is_zero = PETSC_TRUE;
15575ac14e1cSStefano Zampini       }
15585ac14e1cSStefano Zampini     }
15595ac14e1cSStefano Zampini   }
15604cb006feSStefano Zampini   PetscFunctionReturn(0);
15614cb006feSStefano Zampini }
15624cb006feSStefano Zampini 
15634cb006feSStefano Zampini /*@
1564f1580f4eSBarry Smith    PCHYPRESetAlphaPoissonMatrix - Set vector Poisson matrix for `PCHYPRE` of type ams
15654cb006feSStefano Zampini 
1566f1580f4eSBarry Smith    Collective on pc
15674cb006feSStefano Zampini 
15684cb006feSStefano Zampini    Input Parameters:
15694cb006feSStefano Zampini +  pc - the preconditioning context
15704cb006feSStefano Zampini -  A - the matrix
15714cb006feSStefano Zampini 
15724cb006feSStefano Zampini    Level: intermediate
15734cb006feSStefano Zampini 
1574f1580f4eSBarry Smith    Note:
157595452b02SPatrick Sanan     A should be obtained by discretizing the vector valued Poisson problem with linear finite elements
15764cb006feSStefano Zampini 
1577f1580f4eSBarry Smith    Developer Note:
1578f1580f4eSBarry Smith    This automatically converts the matrix to `MATHYPRE` if it is not already of that type
1579f1580f4eSBarry Smith 
1580f1580f4eSBarry Smith    If this is only for  `PCHYPRE` type of ams it should be called `PCHYPREAMSSetAlphaPoissonMatrix()`
1581f1580f4eSBarry Smith 
1582f1580f4eSBarry Smith .seealso: `PCHYPRE`, `PCHYPRESetDiscreteGradient()`, `PCHYPRESetDiscreteCurl()`, `PCHYPRESetBetaPoissonMatrix()`
15834cb006feSStefano Zampini @*/
15849371c9d4SSatish Balay PetscErrorCode PCHYPRESetAlphaPoissonMatrix(PC pc, Mat A) {
15854cb006feSStefano Zampini   PetscFunctionBegin;
15864cb006feSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15874cb006feSStefano Zampini   PetscValidHeaderSpecific(A, MAT_CLASSID, 2);
15884cb006feSStefano Zampini   PetscCheckSameComm(pc, 1, A, 2);
1589cac4c232SBarry Smith   PetscTryMethod(pc, "PCHYPRESetPoissonMatrix_C", (PC, Mat, PetscBool), (pc, A, PETSC_TRUE));
15904cb006feSStefano Zampini   PetscFunctionReturn(0);
15914cb006feSStefano Zampini }
15924cb006feSStefano Zampini 
15934cb006feSStefano Zampini /*@
1594f1580f4eSBarry Smith    PCHYPRESetBetaPoissonMatrix - Set Poisson matrix for `PCHYPRE` of type ams
15954cb006feSStefano Zampini 
1596f1580f4eSBarry Smith    Collective on pc
15974cb006feSStefano Zampini 
15984cb006feSStefano Zampini    Input Parameters:
15994cb006feSStefano Zampini +  pc - the preconditioning context
1600f1580f4eSBarry Smith -  A - the matrix, or NULL to turn it off
16014cb006feSStefano Zampini 
16024cb006feSStefano Zampini    Level: intermediate
16034cb006feSStefano Zampini 
1604f1580f4eSBarry Smith    Note:
160595452b02SPatrick Sanan    A should be obtained by discretizing the Poisson problem with linear finite elements.
16064cb006feSStefano Zampini 
1607f1580f4eSBarry Smith    Developer Note:
1608f1580f4eSBarry Smith    This automatically converts the matrix to `MATHYPRE` if it is not already of that type
1609f1580f4eSBarry Smith 
1610f1580f4eSBarry Smith    If this is only for  `PCHYPRE` type of ams it should be called `PCHYPREAMSPCHYPRESetBetaPoissonMatrix()`
1611f1580f4eSBarry Smith 
1612f1580f4eSBarry Smith .seealso: `PCHYPRE`, `PCHYPRESetDiscreteGradient()`, `PCHYPRESetDiscreteCurl()`, `PCHYPRESetAlphaPoissonMatrix()`
16134cb006feSStefano Zampini @*/
16149371c9d4SSatish Balay PetscErrorCode PCHYPRESetBetaPoissonMatrix(PC pc, Mat A) {
16154cb006feSStefano Zampini   PetscFunctionBegin;
16164cb006feSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16174cb006feSStefano Zampini   if (A) {
16184cb006feSStefano Zampini     PetscValidHeaderSpecific(A, MAT_CLASSID, 2);
16194cb006feSStefano Zampini     PetscCheckSameComm(pc, 1, A, 2);
16204cb006feSStefano Zampini   }
1621cac4c232SBarry Smith   PetscTryMethod(pc, "PCHYPRESetPoissonMatrix_C", (PC, Mat, PetscBool), (pc, A, PETSC_FALSE));
16224cb006feSStefano Zampini   PetscFunctionReturn(0);
16234cb006feSStefano Zampini }
16244cb006feSStefano Zampini 
16259371c9d4SSatish Balay static PetscErrorCode PCHYPRESetEdgeConstantVectors_HYPRE(PC pc, Vec ozz, Vec zoz, Vec zzo) {
16264cb006feSStefano Zampini   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
16274cb006feSStefano Zampini 
16284cb006feSStefano Zampini   PetscFunctionBegin;
16294cb006feSStefano Zampini   /* throw away any vector if already set */
16309566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorDestroy(&jac->constants[0]));
16319566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorDestroy(&jac->constants[1]));
16329566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorDestroy(&jac->constants[2]));
16339566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorCreate(ozz->map, &jac->constants[0]));
16349566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorCopy(ozz, jac->constants[0]));
16359566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorCreate(zoz->map, &jac->constants[1]));
16369566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorCopy(zoz, jac->constants[1]));
16375ac14e1cSStefano Zampini   jac->dim = 2;
16384cb006feSStefano Zampini   if (zzo) {
16399566063dSJacob Faibussowitsch     PetscCall(VecHYPRE_IJVectorCreate(zzo->map, &jac->constants[2]));
16409566063dSJacob Faibussowitsch     PetscCall(VecHYPRE_IJVectorCopy(zzo, jac->constants[2]));
16415ac14e1cSStefano Zampini     jac->dim++;
16424cb006feSStefano Zampini   }
16434cb006feSStefano Zampini   PetscFunctionReturn(0);
16444cb006feSStefano Zampini }
16454cb006feSStefano Zampini 
16464cb006feSStefano Zampini /*@
1647f1580f4eSBarry Smith    PCHYPRESetEdgeConstantVectors - Set the representation of the constant vector fields in the edge element basis for `PCHYPRE` of type ams
16484cb006feSStefano Zampini 
1649f1580f4eSBarry Smith    Collective on pc
16504cb006feSStefano Zampini 
16514cb006feSStefano Zampini    Input Parameters:
16524cb006feSStefano Zampini +  pc - the preconditioning context
16534cb006feSStefano Zampini -  ozz - vector representing (1,0,0) (or (1,0) in 2D)
16544cb006feSStefano Zampini -  zoz - vector representing (0,1,0) (or (0,1) in 2D)
16554cb006feSStefano Zampini -  zzo - vector representing (0,0,1) (use NULL in 2D)
16564cb006feSStefano Zampini 
16574cb006feSStefano Zampini    Level: intermediate
16584cb006feSStefano Zampini 
1659f1580f4eSBarry Smith    Developer Note:
1660f1580f4eSBarry Smith    If this is only for  `PCHYPRE` type of ams it should be called `PCHYPREAMSSetEdgeConstantVectors()`
1661f1580f4eSBarry Smith 
1662f1580f4eSBarry Smith .seealso: `PCHYPRE`, `PCHYPRESetDiscreteGradient()`, `PCHYPRESetDiscreteCurl()`, `PCHYPRESetAlphaPoissonMatrix()`
16634cb006feSStefano Zampini @*/
16649371c9d4SSatish Balay PetscErrorCode PCHYPRESetEdgeConstantVectors(PC pc, Vec ozz, Vec zoz, Vec zzo) {
16654cb006feSStefano Zampini   PetscFunctionBegin;
16664cb006feSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16674cb006feSStefano Zampini   PetscValidHeaderSpecific(ozz, VEC_CLASSID, 2);
16684cb006feSStefano Zampini   PetscValidHeaderSpecific(zoz, VEC_CLASSID, 3);
16694cb006feSStefano Zampini   if (zzo) PetscValidHeaderSpecific(zzo, VEC_CLASSID, 4);
16704cb006feSStefano Zampini   PetscCheckSameComm(pc, 1, ozz, 2);
16714cb006feSStefano Zampini   PetscCheckSameComm(pc, 1, zoz, 3);
16724cb006feSStefano Zampini   if (zzo) PetscCheckSameComm(pc, 1, zzo, 4);
1673cac4c232SBarry Smith   PetscTryMethod(pc, "PCHYPRESetEdgeConstantVectors_C", (PC, Vec, Vec, Vec), (pc, ozz, zoz, zzo));
16744cb006feSStefano Zampini   PetscFunctionReturn(0);
16754cb006feSStefano Zampini }
16764cb006feSStefano Zampini 
16779371c9d4SSatish Balay static PetscErrorCode PCHYPREAMSSetInteriorNodes_HYPRE(PC pc, Vec interior) {
1678be14dc20SKerry Key   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
1679be14dc20SKerry Key 
1680be14dc20SKerry Key   PetscFunctionBegin;
1681be14dc20SKerry Key   PetscCall(VecHYPRE_IJVectorDestroy(&jac->interior));
1682be14dc20SKerry Key   PetscCall(VecHYPRE_IJVectorCreate(interior->map, &jac->interior));
1683be14dc20SKerry Key   PetscCall(VecHYPRE_IJVectorCopy(interior, jac->interior));
1684be14dc20SKerry Key   jac->ams_beta_is_zero_part = PETSC_TRUE;
1685be14dc20SKerry Key   PetscFunctionReturn(0);
1686be14dc20SKerry Key }
1687be14dc20SKerry Key 
1688be14dc20SKerry Key /*@
1689f1580f4eSBarry Smith   PCHYPREAMSSetInteriorNodes - Set the list of interior nodes to a zero-conductivity region for `PCHYPRE` of type ams
1690be14dc20SKerry Key 
1691f1580f4eSBarry Smith    Collective on pc
1692be14dc20SKerry Key 
1693be14dc20SKerry Key    Input Parameters:
1694be14dc20SKerry Key +  pc - the preconditioning context
1695be14dc20SKerry Key -  interior - vector. node is interior if its entry in the array is 1.0.
1696be14dc20SKerry Key 
1697be14dc20SKerry Key    Level: intermediate
1698be14dc20SKerry Key 
1699be14dc20SKerry Key    Note:
1700f1580f4eSBarry Smith    This calls `HYPRE_AMSSetInteriorNodes()`
1701f1580f4eSBarry Smith 
1702f1580f4eSBarry Smith    Developer Note:
1703f1580f4eSBarry Smith    If this is only for  `PCHYPRE` type of ams it should be called `PCHYPREAMSSetInteriorNodes()`
1704f1580f4eSBarry Smith 
1705f1580f4eSBarry Smith .seealso: `PCHYPRE`, `PCHYPRESetDiscreteGradient()`, `PCHYPRESetDiscreteCurl()`, `PCHYPRESetAlphaPoissonMatrix()`
1706be14dc20SKerry Key @*/
17079371c9d4SSatish Balay PetscErrorCode PCHYPREAMSSetInteriorNodes(PC pc, Vec interior) {
1708be14dc20SKerry Key   PetscFunctionBegin;
1709be14dc20SKerry Key   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1710be14dc20SKerry Key   PetscValidHeaderSpecific(interior, VEC_CLASSID, 2);
1711be14dc20SKerry Key   PetscCheckSameComm(pc, 1, interior, 2);
1712be14dc20SKerry Key   PetscTryMethod(pc, "PCHYPREAMSSetInteriorNodes_C", (PC, Vec), (pc, interior));
1713be14dc20SKerry Key   PetscFunctionReturn(0);
1714be14dc20SKerry Key }
1715be14dc20SKerry Key 
17169371c9d4SSatish Balay static PetscErrorCode PCSetCoordinates_HYPRE(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords) {
17174cb006feSStefano Zampini   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
17184cb006feSStefano Zampini   Vec       tv;
17194cb006feSStefano Zampini   PetscInt  i;
17204cb006feSStefano Zampini 
17214cb006feSStefano Zampini   PetscFunctionBegin;
17224cb006feSStefano Zampini   /* throw away any coordinate vector if already set */
17239566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorDestroy(&jac->coords[0]));
17249566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorDestroy(&jac->coords[1]));
17259566063dSJacob Faibussowitsch   PetscCall(VecHYPRE_IJVectorDestroy(&jac->coords[2]));
17265ac14e1cSStefano Zampini   jac->dim = dim;
17275ac14e1cSStefano Zampini 
17284cb006feSStefano Zampini   /* compute IJ vector for coordinates */
17299566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)pc), &tv));
17309566063dSJacob Faibussowitsch   PetscCall(VecSetType(tv, VECSTANDARD));
17319566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(tv, nloc, PETSC_DECIDE));
17324cb006feSStefano Zampini   for (i = 0; i < dim; i++) {
17334cb006feSStefano Zampini     PetscScalar *array;
17344cb006feSStefano Zampini     PetscInt     j;
17354cb006feSStefano Zampini 
17369566063dSJacob Faibussowitsch     PetscCall(VecHYPRE_IJVectorCreate(tv->map, &jac->coords[i]));
17379566063dSJacob Faibussowitsch     PetscCall(VecGetArrayWrite(tv, &array));
17386ea7df73SStefano Zampini     for (j = 0; j < nloc; j++) array[j] = coords[j * dim + i];
17399566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayWrite(tv, &array));
17409566063dSJacob Faibussowitsch     PetscCall(VecHYPRE_IJVectorCopy(tv, jac->coords[i]));
17414cb006feSStefano Zampini   }
17429566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tv));
17434cb006feSStefano Zampini   PetscFunctionReturn(0);
17444cb006feSStefano Zampini }
17454cb006feSStefano Zampini 
17469371c9d4SSatish Balay static PetscErrorCode PCHYPREGetType_HYPRE(PC pc, const char *name[]) {
174716d9e3a6SLisandro Dalcin   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
174816d9e3a6SLisandro Dalcin 
174916d9e3a6SLisandro Dalcin   PetscFunctionBegin;
175016d9e3a6SLisandro Dalcin   *name = jac->hypre_type;
175116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
175216d9e3a6SLisandro Dalcin }
175316d9e3a6SLisandro Dalcin 
17549371c9d4SSatish Balay static PetscErrorCode PCHYPRESetType_HYPRE(PC pc, const char name[]) {
175516d9e3a6SLisandro Dalcin   PC_HYPRE *jac = (PC_HYPRE *)pc->data;
1756ace3abfcSBarry Smith   PetscBool flag;
175716d9e3a6SLisandro Dalcin 
175816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
175916d9e3a6SLisandro Dalcin   if (jac->hypre_type) {
17609566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(jac->hypre_type, name, &flag));
17615f80ce2aSJacob Faibussowitsch     PetscCheck(flag, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Cannot reset the HYPRE preconditioner type once it has been set");
176216d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
176316d9e3a6SLisandro Dalcin   } else {
17649566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(name, &jac->hypre_type));
176516d9e3a6SLisandro Dalcin   }
176616d9e3a6SLisandro Dalcin 
176716d9e3a6SLisandro Dalcin   jac->maxiter         = PETSC_DEFAULT;
176816d9e3a6SLisandro Dalcin   jac->tol             = PETSC_DEFAULT;
176916d9e3a6SLisandro Dalcin   jac->printstatistics = PetscLogPrintInfo;
177016d9e3a6SLisandro Dalcin 
17719566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("pilut", jac->hypre_type, &flag));
177216d9e3a6SLisandro Dalcin   if (flag) {
17739566063dSJacob Faibussowitsch     PetscCall(PetscCommGetComm(PetscObjectComm((PetscObject)pc), &jac->comm_hypre));
1774792fecdfSBarry Smith     PetscCallExternal(HYPRE_ParCSRPilutCreate, jac->comm_hypre, &jac->hsolver);
177516d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
177616d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Pilut;
177716d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParCSRPilutDestroy;
177816d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParCSRPilutSetup;
177916d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParCSRPilutSolve;
178016d9e3a6SLisandro Dalcin     jac->factorrowsize      = PETSC_DEFAULT;
178116d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
178216d9e3a6SLisandro Dalcin   }
17839566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("euclid", jac->hypre_type, &flag));
1784db966c6cSHong Zhang   if (flag) {
17854e3c431bSBarry Smith #if defined(PETSC_USE_64BIT_INDICES)
17864e3c431bSBarry Smith     SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Hypre Euclid does not support 64 bit indices");
17878bf83915SBarry Smith #endif
17889566063dSJacob Faibussowitsch     PetscCall(PetscCommGetComm(PetscObjectComm((PetscObject)pc), &jac->comm_hypre));
1789792fecdfSBarry Smith     PetscCallExternal(HYPRE_EuclidCreate, jac->comm_hypre, &jac->hsolver);
1790db966c6cSHong Zhang     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
1791db966c6cSHong Zhang     pc->ops->view           = PCView_HYPRE_Euclid;
1792db966c6cSHong Zhang     jac->destroy            = HYPRE_EuclidDestroy;
1793db966c6cSHong Zhang     jac->setup              = HYPRE_EuclidSetup;
1794db966c6cSHong Zhang     jac->solve              = HYPRE_EuclidSolve;
1795db966c6cSHong Zhang     jac->factorrowsize      = PETSC_DEFAULT;
1796db966c6cSHong Zhang     jac->eu_level           = PETSC_DEFAULT; /* default */
1797db966c6cSHong Zhang     PetscFunctionReturn(0);
1798db966c6cSHong Zhang   }
17999566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("parasails", jac->hypre_type, &flag));
180016d9e3a6SLisandro Dalcin   if (flag) {
18019566063dSJacob Faibussowitsch     PetscCall(PetscCommGetComm(PetscObjectComm((PetscObject)pc), &jac->comm_hypre));
1802792fecdfSBarry Smith     PetscCallExternal(HYPRE_ParaSailsCreate, jac->comm_hypre, &jac->hsolver);
180316d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
180416d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_ParaSails;
180516d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParaSailsDestroy;
180616d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParaSailsSetup;
180716d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParaSailsSolve;
180816d9e3a6SLisandro Dalcin     /* initialize */
180916d9e3a6SLisandro Dalcin     jac->nlevels            = 1;
18108966356dSPierre Jolivet     jac->threshold          = .1;
181116d9e3a6SLisandro Dalcin     jac->filter             = .1;
181216d9e3a6SLisandro Dalcin     jac->loadbal            = 0;
18132fa5cd67SKarl Rupp     if (PetscLogPrintInfo) jac->logging = (int)PETSC_TRUE;
18142fa5cd67SKarl Rupp     else jac->logging = (int)PETSC_FALSE;
18152fa5cd67SKarl Rupp 
181616d9e3a6SLisandro Dalcin     jac->ruse = (int)PETSC_FALSE;
181716d9e3a6SLisandro Dalcin     jac->symt = 0;
1818792fecdfSBarry Smith     PetscCallExternal(HYPRE_ParaSailsSetParams, jac->hsolver, jac->threshold, jac->nlevels);
1819792fecdfSBarry Smith     PetscCallExternal(HYPRE_ParaSailsSetFilter, jac->hsolver, jac->filter);
1820792fecdfSBarry Smith     PetscCallExternal(HYPRE_ParaSailsSetLoadbal, jac->hsolver, jac->loadbal);
1821792fecdfSBarry Smith     PetscCallExternal(HYPRE_ParaSailsSetLogging, jac->hsolver, jac->logging);
1822792fecdfSBarry Smith     PetscCallExternal(HYPRE_ParaSailsSetReuse, jac->hsolver, jac->ruse);
1823792fecdfSBarry Smith     PetscCallExternal(HYPRE_ParaSailsSetSym, jac->hsolver, jac->symt);
182416d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
182516d9e3a6SLisandro Dalcin   }
18269566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("boomeramg", jac->hypre_type, &flag));
182716d9e3a6SLisandro Dalcin   if (flag) {
1828792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGCreate, &jac->hsolver);
182916d9e3a6SLisandro Dalcin     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
183016d9e3a6SLisandro Dalcin     pc->ops->view            = PCView_HYPRE_BoomerAMG;
183116d9e3a6SLisandro Dalcin     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
183216d9e3a6SLisandro Dalcin     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
18339566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", PCGetInterpolations_BoomerAMG));
18349566063dSJacob Faibussowitsch     PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", PCGetCoarseOperators_BoomerAMG));
183516d9e3a6SLisandro Dalcin     jac->destroy         = HYPRE_BoomerAMGDestroy;
183616d9e3a6SLisandro Dalcin     jac->setup           = HYPRE_BoomerAMGSetup;
183716d9e3a6SLisandro Dalcin     jac->solve           = HYPRE_BoomerAMGSolve;
183816d9e3a6SLisandro Dalcin     jac->applyrichardson = PETSC_FALSE;
183916d9e3a6SLisandro Dalcin     /* these defaults match the hypre defaults */
184016d9e3a6SLisandro Dalcin     jac->cycletype       = 1;
184116d9e3a6SLisandro Dalcin     jac->maxlevels       = 25;
184216d9e3a6SLisandro Dalcin     jac->maxiter         = 1;
18438f87f92bSBarry Smith     jac->tol             = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */
184416d9e3a6SLisandro Dalcin     jac->truncfactor     = 0.0;
184516d9e3a6SLisandro Dalcin     jac->strongthreshold = .25;
184616d9e3a6SLisandro Dalcin     jac->maxrowsum       = .9;
184716d9e3a6SLisandro Dalcin     jac->coarsentype     = 6;
184816d9e3a6SLisandro Dalcin     jac->measuretype     = 0;
18490f1074feSSatish Balay     jac->gridsweeps[0] = jac->gridsweeps[1] = jac->gridsweeps[2] = 1;
18506a251517SEike Mueller     jac->smoothtype                                              = -1; /* Not set by default */
1851b9eb5777SEike Mueller     jac->smoothnumlevels                                         = 25;
18521810e44eSEike Mueller     jac->eu_level                                                = 0;
18531810e44eSEike Mueller     jac->eu_droptolerance                                        = 0;
18541810e44eSEike Mueller     jac->eu_bj                                                   = 0;
1855589dcaf0SStefano Zampini     jac->relaxtype[0] = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a PC - most likely with CG */
18560f1074feSSatish Balay     jac->relaxtype[2]                     = 9; /*G.E. */
185716d9e3a6SLisandro Dalcin     jac->relaxweight                      = 1.0;
185816d9e3a6SLisandro Dalcin     jac->outerrelaxweight                 = 1.0;
185916d9e3a6SLisandro Dalcin     jac->relaxorder                       = 1;
18600f1074feSSatish Balay     jac->interptype                       = 0;
1861589dcaf0SStefano Zampini     jac->Rtype                            = 0;
1862589dcaf0SStefano Zampini     jac->Rstrongthreshold                 = 0.25;
1863589dcaf0SStefano Zampini     jac->Rfilterthreshold                 = 0.0;
1864589dcaf0SStefano Zampini     jac->Adroptype                        = -1;
1865589dcaf0SStefano Zampini     jac->Adroptol                         = 0.0;
18660f1074feSSatish Balay     jac->agg_nl                           = 0;
18676ea7df73SStefano Zampini     jac->agg_interptype                   = 4;
18680f1074feSSatish Balay     jac->pmax                             = 0;
18690f1074feSSatish Balay     jac->truncfactor                      = 0.0;
18700f1074feSSatish Balay     jac->agg_num_paths                    = 1;
1871589dcaf0SStefano Zampini     jac->maxc                             = 9;
1872589dcaf0SStefano Zampini     jac->minc                             = 1;
187322e51d31SStefano Zampini     jac->nodal_coarsening                 = 0;
187422e51d31SStefano Zampini     jac->nodal_coarsening_diag            = 0;
187522e51d31SStefano Zampini     jac->vec_interp_variant               = 0;
187622e51d31SStefano Zampini     jac->vec_interp_qmax                  = 0;
187722e51d31SStefano Zampini     jac->vec_interp_smooth                = PETSC_FALSE;
187822e51d31SStefano Zampini     jac->interp_refine                    = 0;
18798f87f92bSBarry Smith     jac->nodal_relax                      = PETSC_FALSE;
18808f87f92bSBarry Smith     jac->nodal_relax_levels               = 1;
18816ea7df73SStefano Zampini     jac->rap2                             = 0;
18826ea7df73SStefano Zampini 
18836ea7df73SStefano Zampini     /* GPU defaults
18846ea7df73SStefano Zampini          from https://hypre.readthedocs.io/en/latest/solvers-boomeramg.html#gpu-supported-options
18856ea7df73SStefano Zampini          and /src/parcsr_ls/par_amg.c */
18866ea7df73SStefano Zampini #if defined(PETSC_HAVE_HYPRE_DEVICE)
18876ea7df73SStefano Zampini     jac->keeptranspose  = PETSC_TRUE;
18886ea7df73SStefano Zampini     jac->mod_rap2       = 1;
18896ea7df73SStefano Zampini     jac->coarsentype    = 8;
18906ea7df73SStefano Zampini     jac->relaxorder     = 0;
18916ea7df73SStefano Zampini     jac->interptype     = 6;
18926ea7df73SStefano Zampini     jac->relaxtype[0]   = 18;
18936ea7df73SStefano Zampini     jac->relaxtype[1]   = 18;
18946ea7df73SStefano Zampini     jac->agg_interptype = 7;
18956ea7df73SStefano Zampini #else
18966ea7df73SStefano Zampini     jac->keeptranspose = PETSC_FALSE;
18976ea7df73SStefano Zampini     jac->mod_rap2      = 0;
18986ea7df73SStefano Zampini #endif
1899792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetCycleType, jac->hsolver, jac->cycletype);
1900792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetMaxLevels, jac->hsolver, jac->maxlevels);
1901792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetMaxIter, jac->hsolver, jac->maxiter);
1902792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetTol, jac->hsolver, jac->tol);
1903792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetTruncFactor, jac->hsolver, jac->truncfactor);
1904792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetStrongThreshold, jac->hsolver, jac->strongthreshold);
1905792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetMaxRowSum, jac->hsolver, jac->maxrowsum);
1906792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetCoarsenType, jac->hsolver, jac->coarsentype);
1907792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetMeasureType, jac->hsolver, jac->measuretype);
1908792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetRelaxOrder, jac->hsolver, jac->relaxorder);
1909792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetInterpType, jac->hsolver, jac->interptype);
1910792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetAggNumLevels, jac->hsolver, jac->agg_nl);
1911792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetAggInterpType, jac->hsolver, jac->agg_interptype);
1912792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetPMaxElmts, jac->hsolver, jac->pmax);
1913792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetNumPaths, jac->hsolver, jac->agg_num_paths);
1914792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetRelaxType, jac->hsolver, jac->relaxtype[0]);  /* defaults coarse to 9 */
1915792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetNumSweeps, jac->hsolver, jac->gridsweeps[0]); /* defaults coarse to 1 */
1916792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetMaxCoarseSize, jac->hsolver, jac->maxc);
1917792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetMinCoarseSize, jac->hsolver, jac->minc);
19186ea7df73SStefano Zampini     /* GPU */
19196ea7df73SStefano Zampini #if PETSC_PKG_HYPRE_VERSION_GE(2, 18, 0)
1920792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetKeepTranspose, jac->hsolver, jac->keeptranspose ? 1 : 0);
1921792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetRAP2, jac->hsolver, jac->rap2);
1922792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetModuleRAP2, jac->hsolver, jac->mod_rap2);
19236ea7df73SStefano Zampini #endif
19246ea7df73SStefano Zampini 
1925589dcaf0SStefano Zampini     /* AIR */
19266ea7df73SStefano Zampini #if PETSC_PKG_HYPRE_VERSION_GE(2, 18, 0)
1927792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetRestriction, jac->hsolver, jac->Rtype);
1928792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetStrongThresholdR, jac->hsolver, jac->Rstrongthreshold);
1929792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetFilterThresholdR, jac->hsolver, jac->Rfilterthreshold);
1930792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetADropTol, jac->hsolver, jac->Adroptol);
1931792fecdfSBarry Smith     PetscCallExternal(HYPRE_BoomerAMGSetADropType, jac->hsolver, jac->Adroptype);
19326ea7df73SStefano Zampini #endif
193316d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
193416d9e3a6SLisandro Dalcin   }
19359566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("ams", jac->hypre_type, &flag));
19364cb006feSStefano Zampini   if (flag) {
19379566063dSJacob Faibussowitsch     PetscCall(HYPRE_AMSCreate(&jac->hsolver));
19384cb006feSStefano Zampini     pc->ops->setfromoptions   = PCSetFromOptions_HYPRE_AMS;
19394cb006feSStefano Zampini     pc->ops->view             = PCView_HYPRE_AMS;
19404cb006feSStefano Zampini     jac->destroy              = HYPRE_AMSDestroy;
19414cb006feSStefano Zampini     jac->setup                = HYPRE_AMSSetup;
19424cb006feSStefano Zampini     jac->solve                = HYPRE_AMSSolve;
19434cb006feSStefano Zampini     jac->coords[0]            = NULL;
19444cb006feSStefano Zampini     jac->coords[1]            = NULL;
19454cb006feSStefano Zampini     jac->coords[2]            = NULL;
1946be14dc20SKerry Key     jac->interior             = NULL;
19474cb006feSStefano Zampini     /* solver parameters: these are borrowed from mfem package, and they are not the default values from HYPRE AMS */
1948863406b8SStefano Zampini     jac->as_print             = 0;
1949863406b8SStefano Zampini     jac->as_max_iter          = 1;  /* used as a preconditioner */
1950863406b8SStefano Zampini     jac->as_tol               = 0.; /* used as a preconditioner */
19514cb006feSStefano Zampini     jac->ams_cycle_type       = 13;
19524cb006feSStefano Zampini     /* Smoothing options */
1953863406b8SStefano Zampini     jac->as_relax_type        = 2;
1954863406b8SStefano Zampini     jac->as_relax_times       = 1;
1955863406b8SStefano Zampini     jac->as_relax_weight      = 1.0;
1956863406b8SStefano Zampini     jac->as_omega             = 1.0;
19574cb006feSStefano Zampini     /* Vector valued Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
1958863406b8SStefano Zampini     jac->as_amg_alpha_opts[0] = 10;
1959863406b8SStefano Zampini     jac->as_amg_alpha_opts[1] = 1;
19600bdd8552SBarry Smith     jac->as_amg_alpha_opts[2] = 6;
1961863406b8SStefano Zampini     jac->as_amg_alpha_opts[3] = 6;
1962863406b8SStefano Zampini     jac->as_amg_alpha_opts[4] = 4;
1963863406b8SStefano Zampini     jac->as_amg_alpha_theta   = 0.25;
19644cb006feSStefano Zampini     /* Scalar Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
1965863406b8SStefano Zampini     jac->as_amg_beta_opts[0]  = 10;
1966863406b8SStefano Zampini     jac->as_amg_beta_opts[1]  = 1;
19670bdd8552SBarry Smith     jac->as_amg_beta_opts[2]  = 6;
1968863406b8SStefano Zampini     jac->as_amg_beta_opts[3]  = 6;
1969863406b8SStefano Zampini     jac->as_amg_beta_opts[4]  = 4;
1970863406b8SStefano Zampini     jac->as_amg_beta_theta    = 0.25;
1971792fecdfSBarry Smith     PetscCallExternal(HYPRE_AMSSetPrintLevel, jac->hsolver, jac->as_print);
1972792fecdfSBarry Smith     PetscCallExternal(HYPRE_AMSSetMaxIter, jac->hsolver, jac->as_max_iter);
1973792fecdfSBarry Smith     PetscCallExternal(HYPRE_AMSSetCycleType, jac->hsolver, jac->ams_cycle_type);
1974792fecdfSBarry Smith     PetscCallExternal(HYPRE_AMSSetTol, jac->hsolver, jac->as_tol);
19759371c9d4SSatish Balay     PetscCallExternal(HYPRE_AMSSetSmoothingOptions, jac->hsolver, jac->as_relax_type, jac->as_relax_times, jac->as_relax_weight, jac->as_omega);
1976792fecdfSBarry Smith     PetscCallExternal(HYPRE_AMSSetAlphaAMGOptions, jac->hsolver, jac->as_amg_alpha_opts[0], /* AMG coarsen type */
1977863406b8SStefano Zampini                       jac->as_amg_alpha_opts[1],                                            /* AMG agg_levels */
1978863406b8SStefano Zampini                       jac->as_amg_alpha_opts[2],                                            /* AMG relax_type */
19799371c9d4SSatish Balay                       jac->as_amg_alpha_theta, jac->as_amg_alpha_opts[3],                   /* AMG interp_type */
1980a74df02fSJacob Faibussowitsch                       jac->as_amg_alpha_opts[4]);                                           /* AMG Pmax */
1981792fecdfSBarry Smith     PetscCallExternal(HYPRE_AMSSetBetaAMGOptions, jac->hsolver, jac->as_amg_beta_opts[0],   /* AMG coarsen type */
1982863406b8SStefano Zampini                       jac->as_amg_beta_opts[1],                                             /* AMG agg_levels */
1983863406b8SStefano Zampini                       jac->as_amg_beta_opts[2],                                             /* AMG relax_type */
19849371c9d4SSatish Balay                       jac->as_amg_beta_theta, jac->as_amg_beta_opts[3],                     /* AMG interp_type */
1985a74df02fSJacob Faibussowitsch                       jac->as_amg_beta_opts[4]);                                            /* AMG Pmax */
198623df4f25SStefano Zampini     /* Zero conductivity */
198723df4f25SStefano Zampini     jac->ams_beta_is_zero      = PETSC_FALSE;
198823df4f25SStefano Zampini     jac->ams_beta_is_zero_part = PETSC_FALSE;
19894cb006feSStefano Zampini     PetscFunctionReturn(0);
19904cb006feSStefano Zampini   }
19919566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("ads", jac->hypre_type, &flag));
1992863406b8SStefano Zampini   if (flag) {
19939566063dSJacob Faibussowitsch     PetscCall(HYPRE_ADSCreate(&jac->hsolver));
1994863406b8SStefano Zampini     pc->ops->setfromoptions   = PCSetFromOptions_HYPRE_ADS;
1995863406b8SStefano Zampini     pc->ops->view             = PCView_HYPRE_ADS;
1996863406b8SStefano Zampini     jac->destroy              = HYPRE_ADSDestroy;
1997863406b8SStefano Zampini     jac->setup                = HYPRE_ADSSetup;
1998863406b8SStefano Zampini     jac->solve                = HYPRE_ADSSolve;
1999863406b8SStefano Zampini     jac->coords[0]            = NULL;
2000863406b8SStefano Zampini     jac->coords[1]            = NULL;
2001863406b8SStefano Zampini     jac->coords[2]            = NULL;
2002863406b8SStefano Zampini     /* solver parameters: these are borrowed from mfem package, and they are not the default values from HYPRE ADS */
2003863406b8SStefano Zampini     jac->as_print             = 0;
2004863406b8SStefano Zampini     jac->as_max_iter          = 1;  /* used as a preconditioner */
2005863406b8SStefano Zampini     jac->as_tol               = 0.; /* used as a preconditioner */
2006863406b8SStefano Zampini     jac->ads_cycle_type       = 13;
2007863406b8SStefano Zampini     /* Smoothing options */
2008863406b8SStefano Zampini     jac->as_relax_type        = 2;
2009863406b8SStefano Zampini     jac->as_relax_times       = 1;
2010863406b8SStefano Zampini     jac->as_relax_weight      = 1.0;
2011863406b8SStefano Zampini     jac->as_omega             = 1.0;
2012863406b8SStefano Zampini     /* AMS solver parameters: cycle_type, coarsen type, agg_levels, relax_type, interp_type, Pmax */
2013863406b8SStefano Zampini     jac->ams_cycle_type       = 14;
2014863406b8SStefano Zampini     jac->as_amg_alpha_opts[0] = 10;
2015863406b8SStefano Zampini     jac->as_amg_alpha_opts[1] = 1;
2016863406b8SStefano Zampini     jac->as_amg_alpha_opts[2] = 6;
2017863406b8SStefano Zampini     jac->as_amg_alpha_opts[3] = 6;
2018863406b8SStefano Zampini     jac->as_amg_alpha_opts[4] = 4;
2019863406b8SStefano Zampini     jac->as_amg_alpha_theta   = 0.25;
2020863406b8SStefano Zampini     /* Vector Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
2021863406b8SStefano Zampini     jac->as_amg_beta_opts[0]  = 10;
2022863406b8SStefano Zampini     jac->as_amg_beta_opts[1]  = 1;
2023863406b8SStefano Zampini     jac->as_amg_beta_opts[2]  = 6;
2024863406b8SStefano Zampini     jac->as_amg_beta_opts[3]  = 6;
2025863406b8SStefano Zampini     jac->as_amg_beta_opts[4]  = 4;
2026863406b8SStefano Zampini     jac->as_amg_beta_theta    = 0.25;
2027792fecdfSBarry Smith     PetscCallExternal(HYPRE_ADSSetPrintLevel, jac->hsolver, jac->as_print);
2028792fecdfSBarry Smith     PetscCallExternal(HYPRE_ADSSetMaxIter, jac->hsolver, jac->as_max_iter);
2029792fecdfSBarry Smith     PetscCallExternal(HYPRE_ADSSetCycleType, jac->hsolver, jac->ams_cycle_type);
2030792fecdfSBarry Smith     PetscCallExternal(HYPRE_ADSSetTol, jac->hsolver, jac->as_tol);
20319371c9d4SSatish Balay     PetscCallExternal(HYPRE_ADSSetSmoothingOptions, jac->hsolver, jac->as_relax_type, jac->as_relax_times, jac->as_relax_weight, jac->as_omega);
2032792fecdfSBarry Smith     PetscCallExternal(HYPRE_ADSSetAMSOptions, jac->hsolver, jac->ams_cycle_type,      /* AMG coarsen type */
2033863406b8SStefano Zampini                       jac->as_amg_alpha_opts[0],                                      /* AMG coarsen type */
2034863406b8SStefano Zampini                       jac->as_amg_alpha_opts[1],                                      /* AMG agg_levels */
2035863406b8SStefano Zampini                       jac->as_amg_alpha_opts[2],                                      /* AMG relax_type */
20369371c9d4SSatish Balay                       jac->as_amg_alpha_theta, jac->as_amg_alpha_opts[3],             /* AMG interp_type */
2037a74df02fSJacob Faibussowitsch                       jac->as_amg_alpha_opts[4]);                                     /* AMG Pmax */
2038792fecdfSBarry Smith     PetscCallExternal(HYPRE_ADSSetAMGOptions, jac->hsolver, jac->as_amg_beta_opts[0], /* AMG coarsen type */
2039863406b8SStefano Zampini                       jac->as_amg_beta_opts[1],                                       /* AMG agg_levels */
2040863406b8SStefano Zampini                       jac->as_amg_beta_opts[2],                                       /* AMG relax_type */
20419371c9d4SSatish Balay                       jac->as_amg_beta_theta, jac->as_amg_beta_opts[3],               /* AMG interp_type */
2042a74df02fSJacob Faibussowitsch                       jac->as_amg_beta_opts[4]);                                      /* AMG Pmax */
2043863406b8SStefano Zampini     PetscFunctionReturn(0);
2044863406b8SStefano Zampini   }
20459566063dSJacob Faibussowitsch   PetscCall(PetscFree(jac->hypre_type));
20462fa5cd67SKarl Rupp 
20470298fd71SBarry Smith   jac->hypre_type = NULL;
204898921bdaSJacob Faibussowitsch   SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown HYPRE preconditioner %s; Choices are euclid, pilut, parasails, boomeramg, ams", name);
204916d9e3a6SLisandro Dalcin }
205016d9e3a6SLisandro Dalcin 
205116d9e3a6SLisandro Dalcin /*
205216d9e3a6SLisandro Dalcin     It only gets here if the HYPRE type has not been set before the call to
205316d9e3a6SLisandro Dalcin    ...SetFromOptions() which actually is most of the time
205416d9e3a6SLisandro Dalcin */
20559371c9d4SSatish Balay PetscErrorCode PCSetFromOptions_HYPRE(PC pc, PetscOptionItems *PetscOptionsObject) {
20564ddd07fcSJed Brown   PetscInt    indx;
2057db966c6cSHong Zhang   const char *type[] = {"euclid", "pilut", "parasails", "boomeramg", "ams", "ads"};
2058ace3abfcSBarry Smith   PetscBool   flg;
205916d9e3a6SLisandro Dalcin 
206016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
2061d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "HYPRE preconditioner options");
2062dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-pc_hypre_type", "HYPRE preconditioner type", "PCHYPRESetType", type, PETSC_STATIC_ARRAY_LENGTH(type), "boomeramg", &indx, &flg));
206316d9e3a6SLisandro Dalcin   if (flg) {
20649566063dSJacob Faibussowitsch     PetscCall(PCHYPRESetType_HYPRE(pc, type[indx]));
206502a17cd4SBarry Smith   } else {
20669566063dSJacob Faibussowitsch     PetscCall(PCHYPRESetType_HYPRE(pc, "boomeramg"));
206716d9e3a6SLisandro Dalcin   }
2068dbbe0bcdSBarry Smith   PetscTryTypeMethod(pc, setfromoptions, PetscOptionsObject);
2069d0609cedSBarry Smith   PetscOptionsHeadEnd();
207016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
207116d9e3a6SLisandro Dalcin }
207216d9e3a6SLisandro Dalcin 
207316d9e3a6SLisandro Dalcin /*@C
207416d9e3a6SLisandro Dalcin      PCHYPRESetType - Sets which hypre preconditioner you wish to use
207516d9e3a6SLisandro Dalcin 
207616d9e3a6SLisandro Dalcin    Input Parameters:
207716d9e3a6SLisandro Dalcin +     pc - the preconditioner context
2078db966c6cSHong Zhang -     name - either  euclid, pilut, parasails, boomeramg, ams, ads
207916d9e3a6SLisandro Dalcin 
2080f1580f4eSBarry Smith    Options Database Key:
2081db966c6cSHong Zhang    -pc_hypre_type - One of euclid, pilut, parasails, boomeramg, ams, ads
208216d9e3a6SLisandro Dalcin 
208316d9e3a6SLisandro Dalcin    Level: intermediate
208416d9e3a6SLisandro Dalcin 
2085f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCHYPRE`
208616d9e3a6SLisandro Dalcin @*/
20879371c9d4SSatish Balay PetscErrorCode PCHYPRESetType(PC pc, const char name[]) {
208816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
20890700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
209016d9e3a6SLisandro Dalcin   PetscValidCharPointer(name, 2);
2091cac4c232SBarry Smith   PetscTryMethod(pc, "PCHYPRESetType_C", (PC, const char[]), (pc, name));
209216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
209316d9e3a6SLisandro Dalcin }
209416d9e3a6SLisandro Dalcin 
209516d9e3a6SLisandro Dalcin /*@C
209616d9e3a6SLisandro Dalcin      PCHYPREGetType - Gets which hypre preconditioner you are using
209716d9e3a6SLisandro Dalcin 
209816d9e3a6SLisandro Dalcin    Input Parameter:
209916d9e3a6SLisandro Dalcin .     pc - the preconditioner context
210016d9e3a6SLisandro Dalcin 
210116d9e3a6SLisandro Dalcin    Output Parameter:
2102db966c6cSHong Zhang .     name - either  euclid, pilut, parasails, boomeramg, ams, ads
210316d9e3a6SLisandro Dalcin 
210416d9e3a6SLisandro Dalcin    Level: intermediate
210516d9e3a6SLisandro Dalcin 
2106f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCHYPRESetType()`, `PCType`, `PC`, `PCHYPRE`
210716d9e3a6SLisandro Dalcin @*/
21089371c9d4SSatish Balay PetscErrorCode PCHYPREGetType(PC pc, const char *name[]) {
210916d9e3a6SLisandro Dalcin   PetscFunctionBegin;
21100700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
211116d9e3a6SLisandro Dalcin   PetscValidPointer(name, 2);
2112cac4c232SBarry Smith   PetscTryMethod(pc, "PCHYPREGetType_C", (PC, const char *[]), (pc, name));
211316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
211416d9e3a6SLisandro Dalcin }
211516d9e3a6SLisandro Dalcin 
2116db6f9c32SMark Adams /*@C
2117f1580f4eSBarry Smith    PCMGGalerkinSetMatProductAlgorithm - Set type of SpGEMM for hypre to use on GPUs
2118db6f9c32SMark Adams 
2119f1580f4eSBarry Smith    Logically Collective on pc
2120db6f9c32SMark Adams 
2121db6f9c32SMark Adams    Input Parameters:
2122db6f9c32SMark Adams +  pc - the hypre context
2123db6f9c32SMark Adams -  type - one of 'cusparse', 'hypre'
2124db6f9c32SMark Adams 
2125db6f9c32SMark Adams    Options Database Key:
212667b8a455SSatish Balay .  -pc_mg_galerkin_mat_product_algorithm <cusparse,hypre> - Type of SpGEMM to use in hypre
2127db6f9c32SMark Adams 
2128db6f9c32SMark Adams    Level: intermediate
2129db6f9c32SMark Adams 
2130f1580f4eSBarry Smith    Developer Note:
2131f1580f4eSBarry Smith    How the name starts with `PCMG`, should it not be `PCHYPREBoomerAMG`?
2132db6f9c32SMark Adams 
2133f1580f4eSBarry Smith .seealso: `PCHYPRE`, `PCMGGalerkinGetMatProductAlgorithm()`
2134db6f9c32SMark Adams @*/
21359371c9d4SSatish Balay PetscErrorCode PCMGGalerkinSetMatProductAlgorithm(PC pc, const char name[]) {
2136db6f9c32SMark Adams   PetscFunctionBegin;
2137db6f9c32SMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2138cac4c232SBarry Smith   PetscTryMethod(pc, "PCMGGalerkinSetMatProductAlgorithm_C", (PC, const char[]), (pc, name));
2139db6f9c32SMark Adams   PetscFunctionReturn(0);
2140db6f9c32SMark Adams }
2141db6f9c32SMark Adams 
2142db6f9c32SMark Adams /*@C
2143f1580f4eSBarry Smith    PCMGGalerkinGetMatProductAlgorithm - Get type of SpGEMM for hypre to use on GPUs
2144db6f9c32SMark Adams 
2145db6f9c32SMark Adams    Not Collective
2146db6f9c32SMark Adams 
2147db6f9c32SMark Adams    Input Parameter:
2148db6f9c32SMark Adams .  pc - the multigrid context
2149db6f9c32SMark Adams 
2150db6f9c32SMark Adams    Output Parameter:
2151db6f9c32SMark Adams .  name - one of 'cusparse', 'hypre'
2152db6f9c32SMark Adams 
2153db6f9c32SMark Adams    Level: intermediate
2154db6f9c32SMark Adams 
2155f1580f4eSBarry Smith .seealso: `PCHYPRE`, ``PCMGGalerkinSetMatProductAlgorithm()`
2156db6f9c32SMark Adams @*/
21579371c9d4SSatish Balay PetscErrorCode PCMGGalerkinGetMatProductAlgorithm(PC pc, const char *name[]) {
2158db6f9c32SMark Adams   PetscFunctionBegin;
2159db6f9c32SMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2160cac4c232SBarry Smith   PetscTryMethod(pc, "PCMGGalerkinGetMatProductAlgorithm_C", (PC, const char *[]), (pc, name));
2161db6f9c32SMark Adams   PetscFunctionReturn(0);
2162db6f9c32SMark Adams }
2163db6f9c32SMark Adams 
216416d9e3a6SLisandro Dalcin /*MC
2165f1580f4eSBarry Smith      PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre as PETSc `PC`
216616d9e3a6SLisandro Dalcin 
216716d9e3a6SLisandro Dalcin    Options Database Keys:
2168db966c6cSHong Zhang +   -pc_hypre_type - One of euclid, pilut, parasails, boomeramg, ams, ads
2169f1580f4eSBarry Smith .   -pc_hypre_boomeramg_nodal_coarsen <n> - where n is from 1 to 6 (see `HYPRE_BOOMERAMGSetNodal()`)
2170f1580f4eSBarry Smith .   -pc_hypre_boomeramg_vec_interp_variant <v> - where v is from 1 to 3 (see `HYPRE_BoomerAMGSetInterpVecVariant()`)
2171ead8c081SBarry Smith -   Many others, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX preconditioner
217216d9e3a6SLisandro Dalcin 
217316d9e3a6SLisandro Dalcin    Level: intermediate
217416d9e3a6SLisandro Dalcin 
217595452b02SPatrick Sanan    Notes:
2176f1580f4eSBarry Smith     Apart from pc_hypre_type (for which there is `PCHYPRESetType()`),
217716d9e3a6SLisandro Dalcin           the many hypre options can ONLY be set via the options database (e.g. the command line
217816d9e3a6SLisandro Dalcin           or with PetscOptionsSetValue(), there are no functions to set them)
217916d9e3a6SLisandro Dalcin 
2180c231f9e3SBarryFSmith           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_tol refer to the number of iterations
21810f1074feSSatish Balay           (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if
21820f1074feSSatish Balay           -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner
2183c231f9e3SBarryFSmith           (-pc_hypre_boomeramg_tol should be set to 0.0 - the default - to strictly use a fixed number of
21848f87f92bSBarry Smith           iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations
21850f1074feSSatish Balay           and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10
21860f1074feSSatish Balay           then AT MOST twenty V-cycles of boomeramg will be called.
218716d9e3a6SLisandro Dalcin 
21880f1074feSSatish Balay            Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation
21890f1074feSSatish Balay            (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry.
21900f1074feSSatish Balay            Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi.
219116d9e3a6SLisandro Dalcin           If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly
219216d9e3a6SLisandro Dalcin           and use -ksp_max_it to control the number of V-cycles.
219316d9e3a6SLisandro Dalcin           (see the PETSc FAQ.html at the PETSc website under the Documentation tab).
219416d9e3a6SLisandro Dalcin 
219516d9e3a6SLisandro Dalcin           2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option
219616d9e3a6SLisandro Dalcin           -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L.
219716d9e3a6SLisandro Dalcin 
2198f1580f4eSBarry Smith           `MatSetNearNullSpace()` - if you provide a near null space to your matrix it is ignored by hypre UNLESS you also use
2199fdd15c9aSJunchao Zhang           the following two options:
22000b1a5bd9SEric Chamberland 
2201f1580f4eSBarry Smith           See `PCPFMG`, `PCSMG`, and `PCSYSPFMG` for access to hypre's other (nonalgebraic) multigrid solvers
2202f1580f4eSBarry Smith 
2203f1580f4eSBarry Smith           For `PCHYPRE` type of ams or ads auxiliary data must be provided to the preconditioner with `PCHYPRESetDiscreteGradient()`,
2204f1580f4eSBarry Smith           `PCHYPRESetDiscreteCurl()`, `PCHYPRESetInterpolations()`, `PCHYPRESetAlphaPoissonMatrix()`, `PCHYPRESetBetaPoissonMatrix()`, `PCHYPRESetEdgeConstantVectors()`,
2205f1580f4eSBarry Smith           PCHYPREAMSSetInteriorNodes()
2206f1580f4eSBarry Smith 
2207f1580f4eSBarry Smith    PETSc provides its own geometric and algebraic multigrid solvers `PCMG` and `PCGAMG`, also see `PCHMG` which is useful for certain multicomponent problems
22089e5bc791SBarry Smith 
2209ead8c081SBarry Smith    GPU Notes:
2210ead8c081SBarry Smith      To configure hypre BoomerAMG so that it can utilize NVIDIA GPUs run ./configure --download-hypre --with-cuda
2211f1580f4eSBarry Smith      Then pass `VECCUDA` vectors and `MATAIJCUSPARSE` matrices to the solvers and PETSc will automatically utilize hypre's GPU solvers.
2212ead8c081SBarry Smith 
2213ead8c081SBarry Smith      To configure hypre BoomerAMG so that it can utilize AMD GPUs run ./configure --download-hypre --with-hip
2214f1580f4eSBarry Smith      Then pass `VECHIP` vectors to the solvers and PETSc will automatically utilize hypre's GPU solvers.
2215ead8c081SBarry Smith 
2216f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCHYPRESetType()`, `PCPFMG`, `PCGAMG`, `PCSYSPFMG`, `PCSMG`, `PCHYPRESetDiscreteGradient()`,
2217f1580f4eSBarry Smith           `PCHYPRESetDiscreteCurl()`, `PCHYPRESetInterpolations()`, `PCHYPRESetAlphaPoissonMatrix()`, `PCHYPRESetBetaPoissonMatrix()`, `PCHYPRESetEdgeConstantVectors()`,
2218f1580f4eSBarry Smith           PCHYPREAMSSetInteriorNodes()
221916d9e3a6SLisandro Dalcin M*/
222016d9e3a6SLisandro Dalcin 
22219371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_HYPRE(PC pc) {
222216d9e3a6SLisandro Dalcin   PC_HYPRE *jac;
222316d9e3a6SLisandro Dalcin 
222416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
2225*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
22262fa5cd67SKarl Rupp 
222716d9e3a6SLisandro Dalcin   pc->data                = jac;
22288695de01SBarry Smith   pc->ops->reset          = PCReset_HYPRE;
222916d9e3a6SLisandro Dalcin   pc->ops->destroy        = PCDestroy_HYPRE;
223016d9e3a6SLisandro Dalcin   pc->ops->setfromoptions = PCSetFromOptions_HYPRE;
223116d9e3a6SLisandro Dalcin   pc->ops->setup          = PCSetUp_HYPRE;
223216d9e3a6SLisandro Dalcin   pc->ops->apply          = PCApply_HYPRE;
223316d9e3a6SLisandro Dalcin   jac->comm_hypre         = MPI_COMM_NULL;
22349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPRESetType_C", PCHYPRESetType_HYPRE));
22359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPREGetType_C", PCHYPREGetType_HYPRE));
22369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetCoordinates_C", PCSetCoordinates_HYPRE));
22379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPRESetDiscreteGradient_C", PCHYPRESetDiscreteGradient_HYPRE));
22389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPRESetDiscreteCurl_C", PCHYPRESetDiscreteCurl_HYPRE));
22399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPRESetInterpolations_C", PCHYPRESetInterpolations_HYPRE));
22409566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPRESetEdgeConstantVectors_C", PCHYPRESetEdgeConstantVectors_HYPRE));
2241be14dc20SKerry Key   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPREAMSSetInteriorNodes_C", PCHYPREAMSSetInteriorNodes_HYPRE));
22429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHYPRESetPoissonMatrix_C", PCHYPRESetPoissonMatrix_HYPRE));
22439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGalerkinSetMatProductAlgorithm_C", PCMGGalerkinSetMatProductAlgorithm_HYPRE_BoomerAMG));
22449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGalerkinGetMatProductAlgorithm_C", PCMGGalerkinGetMatProductAlgorithm_HYPRE_BoomerAMG));
22456ea7df73SStefano Zampini #if defined(PETSC_HAVE_HYPRE_DEVICE)
22466ea7df73SStefano Zampini #if defined(HYPRE_USING_HIP)
22479566063dSJacob Faibussowitsch   PetscCall(PetscDeviceInitialize(PETSC_DEVICE_HIP));
22486ea7df73SStefano Zampini #endif
22496ea7df73SStefano Zampini #if defined(HYPRE_USING_CUDA)
22509566063dSJacob Faibussowitsch   PetscCall(PetscDeviceInitialize(PETSC_DEVICE_CUDA));
22516ea7df73SStefano Zampini #endif
22526ea7df73SStefano Zampini #endif
225316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
225416d9e3a6SLisandro Dalcin }
2255ebc551c0SBarry Smith 
2256ebc551c0SBarry Smith typedef struct {
225768326731SBarry Smith   MPI_Comm           hcomm; /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */
2258f91d8e95SBarry Smith   HYPRE_StructSolver hsolver;
22599e5bc791SBarry Smith 
22609e5bc791SBarry Smith   /* keep copy of PFMG options used so may view them */
22614ddd07fcSJed Brown   PetscInt  its;
22629e5bc791SBarry Smith   double    tol;
22634ddd07fcSJed Brown   PetscInt  relax_type;
22644ddd07fcSJed Brown   PetscInt  rap_type;
22654ddd07fcSJed Brown   PetscInt  num_pre_relax, num_post_relax;
22664ddd07fcSJed Brown   PetscInt  max_levels;
22670be8cd64Sftrigaux   PetscInt  skip_relax;
22680be8cd64Sftrigaux   PetscBool print_statistics;
2269ebc551c0SBarry Smith } PC_PFMG;
2270ebc551c0SBarry Smith 
22719371c9d4SSatish Balay PetscErrorCode PCDestroy_PFMG(PC pc) {
2272f91d8e95SBarry Smith   PC_PFMG *ex = (PC_PFMG *)pc->data;
2273ebc551c0SBarry Smith 
2274ebc551c0SBarry Smith   PetscFunctionBegin;
2275792fecdfSBarry Smith   if (ex->hsolver) PetscCallExternal(HYPRE_StructPFMGDestroy, ex->hsolver);
22769566063dSJacob Faibussowitsch   PetscCall(PetscCommRestoreComm(PetscObjectComm((PetscObject)pc), &ex->hcomm));
22779566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
2278ebc551c0SBarry Smith   PetscFunctionReturn(0);
2279ebc551c0SBarry Smith }
2280ebc551c0SBarry Smith 
22819e5bc791SBarry Smith static const char *PFMGRelaxType[] = {"Jacobi", "Weighted-Jacobi", "symmetric-Red/Black-Gauss-Seidel", "Red/Black-Gauss-Seidel"};
22829e5bc791SBarry Smith static const char *PFMGRAPType[]   = {"Galerkin", "non-Galerkin"};
22839e5bc791SBarry Smith 
22849371c9d4SSatish Balay PetscErrorCode PCView_PFMG(PC pc, PetscViewer viewer) {
2285ace3abfcSBarry Smith   PetscBool iascii;
2286f91d8e95SBarry Smith   PC_PFMG  *ex = (PC_PFMG *)pc->data;
2287ebc551c0SBarry Smith 
2288ebc551c0SBarry Smith   PetscFunctionBegin;
22899566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
22909e5bc791SBarry Smith   if (iascii) {
22919566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  HYPRE PFMG preconditioning\n"));
229263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    max iterations %" PetscInt_FMT "\n", ex->its));
22939566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    tolerance %g\n", ex->tol));
22949566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    relax type %s\n", PFMGRelaxType[ex->relax_type]));
22959566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    RAP type %s\n", PFMGRAPType[ex->rap_type]));
229663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    number pre-relax %" PetscInt_FMT " post-relax %" PetscInt_FMT "\n", ex->num_pre_relax, ex->num_post_relax));
229763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    max levels %" PetscInt_FMT "\n", ex->max_levels));
22980be8cd64Sftrigaux     PetscCall(PetscViewerASCIIPrintf(viewer, "    skip relax %" PetscInt_FMT "\n", ex->skip_relax));
22999e5bc791SBarry Smith   }
2300ebc551c0SBarry Smith   PetscFunctionReturn(0);
2301ebc551c0SBarry Smith }
2302ebc551c0SBarry Smith 
23039371c9d4SSatish Balay PetscErrorCode PCSetFromOptions_PFMG(PC pc, PetscOptionItems *PetscOptionsObject) {
2304f91d8e95SBarry Smith   PC_PFMG *ex = (PC_PFMG *)pc->data;
2305ebc551c0SBarry Smith 
2306ebc551c0SBarry Smith   PetscFunctionBegin;
2307d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "PFMG options");
23080be8cd64Sftrigaux   PetscCall(PetscOptionsBool("-pc_pfmg_print_statistics", "Print statistics", "HYPRE_StructPFMGSetPrintLevel", ex->print_statistics, &ex->print_statistics, NULL));
23099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_pfmg_its", "Number of iterations of PFMG to use as preconditioner", "HYPRE_StructPFMGSetMaxIter", ex->its, &ex->its, NULL));
2310792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGSetMaxIter, ex->hsolver, ex->its);
23119566063dSJacob Faibussowitsch   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));
2312792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGSetNumPreRelax, ex->hsolver, ex->num_pre_relax);
23139566063dSJacob Faibussowitsch   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));
2314792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGSetNumPostRelax, ex->hsolver, ex->num_post_relax);
23159e5bc791SBarry Smith 
23169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_pfmg_max_levels", "Max Levels for MG hierarchy", "HYPRE_StructPFMGSetMaxLevels", ex->max_levels, &ex->max_levels, NULL));
2317792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGSetMaxLevels, ex->hsolver, ex->max_levels);
23183b46a515SGlenn Hammond 
23199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_pfmg_tol", "Tolerance of PFMG", "HYPRE_StructPFMGSetTol", ex->tol, &ex->tol, NULL));
2320792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGSetTol, ex->hsolver, ex->tol);
2321dd39110bSPierre Jolivet   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));
2322792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGSetRelaxType, ex->hsolver, ex->relax_type);
2323dd39110bSPierre Jolivet   PetscCall(PetscOptionsEList("-pc_pfmg_rap_type", "RAP type", "HYPRE_StructPFMGSetRAPType", PFMGRAPType, PETSC_STATIC_ARRAY_LENGTH(PFMGRAPType), PFMGRAPType[ex->rap_type], &ex->rap_type, NULL));
2324792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGSetRAPType, ex->hsolver, ex->rap_type);
23250be8cd64Sftrigaux   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));
23260be8cd64Sftrigaux   PetscCallExternal(HYPRE_StructPFMGSetSkipRelax, ex->hsolver, ex->skip_relax);
2327d0609cedSBarry Smith   PetscOptionsHeadEnd();
2328ebc551c0SBarry Smith   PetscFunctionReturn(0);
2329ebc551c0SBarry Smith }
2330ebc551c0SBarry Smith 
23319371c9d4SSatish Balay PetscErrorCode PCApply_PFMG(PC pc, Vec x, Vec y) {
2332f91d8e95SBarry Smith   PC_PFMG           *ex = (PC_PFMG *)pc->data;
2333d9ca1df4SBarry Smith   PetscScalar       *yy;
2334d9ca1df4SBarry Smith   const PetscScalar *xx;
23354ddd07fcSJed Brown   PetscInt           ilower[3], iupper[3];
23362cf14000SStefano Zampini   HYPRE_Int          hlower[3], hupper[3];
233768326731SBarry Smith   Mat_HYPREStruct   *mx = (Mat_HYPREStruct *)(pc->pmat->data);
2338f91d8e95SBarry Smith 
2339f91d8e95SBarry Smith   PetscFunctionBegin;
23409566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(hypreCitation, &cite));
23419566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(mx->da, &ilower[0], &ilower[1], &ilower[2], &iupper[0], &iupper[1], &iupper[2]));
23422cf14000SStefano Zampini   /* when HYPRE_MIXEDINT is defined, sizeof(HYPRE_Int) == 32 */
2343f91d8e95SBarry Smith   iupper[0] += ilower[0] - 1;
2344f91d8e95SBarry Smith   iupper[1] += ilower[1] - 1;
2345f91d8e95SBarry Smith   iupper[2] += ilower[2] - 1;
23462cf14000SStefano Zampini   hlower[0] = (HYPRE_Int)ilower[0];
23472cf14000SStefano Zampini   hlower[1] = (HYPRE_Int)ilower[1];
23482cf14000SStefano Zampini   hlower[2] = (HYPRE_Int)ilower[2];
23492cf14000SStefano Zampini   hupper[0] = (HYPRE_Int)iupper[0];
23502cf14000SStefano Zampini   hupper[1] = (HYPRE_Int)iupper[1];
23512cf14000SStefano Zampini   hupper[2] = (HYPRE_Int)iupper[2];
2352f91d8e95SBarry Smith 
2353f91d8e95SBarry Smith   /* copy x values over to hypre */
2354792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructVectorSetConstantValues, mx->hb, 0.0);
23559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xx));
2356792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructVectorSetBoxValues, mx->hb, hlower, hupper, (HYPRE_Complex *)xx);
23579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xx));
2358792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructVectorAssemble, mx->hb);
2359792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGSolve, ex->hsolver, mx->hmat, mx->hb, mx->hx);
2360f91d8e95SBarry Smith 
2361f91d8e95SBarry Smith   /* copy solution values back to PETSc */
23629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &yy));
2363792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructVectorGetBoxValues, mx->hx, hlower, hupper, (HYPRE_Complex *)yy);
23649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &yy));
2365f91d8e95SBarry Smith   PetscFunctionReturn(0);
2366f91d8e95SBarry Smith }
2367f91d8e95SBarry Smith 
23689371c9d4SSatish Balay 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) {
23699e5bc791SBarry Smith   PC_PFMG  *jac = (PC_PFMG *)pc->data;
23702cf14000SStefano Zampini   HYPRE_Int oits;
23719e5bc791SBarry Smith 
23729e5bc791SBarry Smith   PetscFunctionBegin;
23739566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(hypreCitation, &cite));
2374792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGSetMaxIter, jac->hsolver, its * jac->its);
2375792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGSetTol, jac->hsolver, rtol);
23769e5bc791SBarry Smith 
23779566063dSJacob Faibussowitsch   PetscCall(PCApply_PFMG(pc, b, y));
2378792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGGetNumIterations, jac->hsolver, &oits);
23799e5bc791SBarry Smith   *outits = oits;
23809e5bc791SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
23819e5bc791SBarry Smith   else *reason = PCRICHARDSON_CONVERGED_RTOL;
2382792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGSetTol, jac->hsolver, jac->tol);
2383792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGSetMaxIter, jac->hsolver, jac->its);
23849e5bc791SBarry Smith   PetscFunctionReturn(0);
23859e5bc791SBarry Smith }
23869e5bc791SBarry Smith 
23879371c9d4SSatish Balay PetscErrorCode PCSetUp_PFMG(PC pc) {
23883a32d3dbSGlenn Hammond   PC_PFMG         *ex = (PC_PFMG *)pc->data;
23893a32d3dbSGlenn Hammond   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
2390ace3abfcSBarry Smith   PetscBool        flg;
23913a32d3dbSGlenn Hammond 
23923a32d3dbSGlenn Hammond   PetscFunctionBegin;
23939566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATHYPRESTRUCT, &flg));
239428b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Must use MATHYPRESTRUCT with this preconditioner");
23953a32d3dbSGlenn Hammond 
23963a32d3dbSGlenn Hammond   /* create the hypre solver object and set its information */
2397792fecdfSBarry Smith   if (ex->hsolver) PetscCallExternal(HYPRE_StructPFMGDestroy, ex->hsolver);
2398792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGCreate, ex->hcomm, &ex->hsolver);
23990be8cd64Sftrigaux 
24000be8cd64Sftrigaux   // Print Hypre statistics about the solve process
24010be8cd64Sftrigaux   if (ex->print_statistics) PetscCallExternal(HYPRE_StructPFMGSetPrintLevel, ex->hsolver, 3);
24020be8cd64Sftrigaux 
24030be8cd64Sftrigaux   // The hypre options must be repeated here because the StructPFMG was destroyed and recreated
24040be8cd64Sftrigaux   PetscCallExternal(HYPRE_StructPFMGSetMaxIter, ex->hsolver, ex->its);
24050be8cd64Sftrigaux   PetscCallExternal(HYPRE_StructPFMGSetNumPreRelax, ex->hsolver, ex->num_pre_relax);
24060be8cd64Sftrigaux   PetscCallExternal(HYPRE_StructPFMGSetNumPostRelax, ex->hsolver, ex->num_post_relax);
24070be8cd64Sftrigaux   PetscCallExternal(HYPRE_StructPFMGSetMaxLevels, ex->hsolver, ex->max_levels);
24080be8cd64Sftrigaux   PetscCallExternal(HYPRE_StructPFMGSetTol, ex->hsolver, ex->tol);
24090be8cd64Sftrigaux   PetscCallExternal(HYPRE_StructPFMGSetRelaxType, ex->hsolver, ex->relax_type);
24100be8cd64Sftrigaux   PetscCallExternal(HYPRE_StructPFMGSetRAPType, ex->hsolver, ex->rap_type);
24110be8cd64Sftrigaux 
2412792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGSetup, ex->hsolver, mx->hmat, mx->hb, mx->hx);
2413792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGSetZeroGuess, ex->hsolver);
24143a32d3dbSGlenn Hammond   PetscFunctionReturn(0);
24153a32d3dbSGlenn Hammond }
24163a32d3dbSGlenn Hammond 
2417ebc551c0SBarry Smith /*MC
2418ebc551c0SBarry Smith      PCPFMG - the hypre PFMG multigrid solver
2419ebc551c0SBarry Smith 
2420f1580f4eSBarry Smith    Options Database Keys:
242167b8a455SSatish Balay + -pc_pfmg_its <its> - number of iterations of PFMG to use as preconditioner
242267b8a455SSatish Balay . -pc_pfmg_num_pre_relax <steps> - number of smoothing steps before coarse grid solve
242367b8a455SSatish Balay . -pc_pfmg_num_post_relax <steps> - number of smoothing steps after coarse grid solve
242467b8a455SSatish Balay . -pc_pfmg_tol <tol> - tolerance of PFMG
24259e5bc791SBarry Smith . -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
24260be8cd64Sftrigaux . -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin
2427f1580f4eSBarry Smith - -pc_pfmg_skip_relax - skip relaxation on certain grids for isotropic problems. This can greatly improve efficiency by eliminating unnecessary relaxations
2428f1580f4eSBarry Smith                         when the underlying problem is isotropic, one of 0,1
2429f1580f4eSBarry Smith 
2430f1580f4eSBarry Smith    Level: advanced
2431f91d8e95SBarry Smith 
243295452b02SPatrick Sanan    Notes:
243395452b02SPatrick Sanan    This is for CELL-centered descretizations
24349e5bc791SBarry Smith 
2435f1580f4eSBarry Smith    See `PCSYSPFMG` for a version suitable for systems of PDEs, and `PCSMG`
24369e5bc791SBarry Smith 
2437f1580f4eSBarry Smith    See `PCHYPRE` for hypre's BoomerAMG algebraic multigrid solver
2438f1580f4eSBarry Smith 
2439f1580f4eSBarry Smith    This must be used with the `MATHYPRESTRUCT` matrix type.
2440f1580f4eSBarry Smith 
2441f1580f4eSBarry Smith    This provides only some of the functionality of PFMG, it supports only one block per process defined by a PETSc `DMDA`.
2442f1580f4eSBarry Smith 
2443f1580f4eSBarry Smith .seealso: `PCMG`, `MATHYPRESTRUCT`, `PCHYPRE`, `PCGAMG`, `PCSYSPFMG`, `PCSMG`
2444ebc551c0SBarry Smith M*/
2445ebc551c0SBarry Smith 
24469371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_PFMG(PC pc) {
2447ebc551c0SBarry Smith   PC_PFMG *ex;
2448ebc551c0SBarry Smith 
2449ebc551c0SBarry Smith   PetscFunctionBegin;
24509371c9d4SSatish Balay   PetscCall(PetscNew(&ex));
245168326731SBarry Smith   pc->data = ex;
2452ebc551c0SBarry Smith 
24539e5bc791SBarry Smith   ex->its              = 1;
24549e5bc791SBarry Smith   ex->tol              = 1.e-8;
24559e5bc791SBarry Smith   ex->relax_type       = 1;
24569e5bc791SBarry Smith   ex->rap_type         = 0;
24579e5bc791SBarry Smith   ex->num_pre_relax    = 1;
24589e5bc791SBarry Smith   ex->num_post_relax   = 1;
24593b46a515SGlenn Hammond   ex->max_levels       = 0;
24600be8cd64Sftrigaux   ex->skip_relax       = 0;
24610be8cd64Sftrigaux   ex->print_statistics = PETSC_FALSE;
24629e5bc791SBarry Smith 
2463ebc551c0SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_PFMG;
2464ebc551c0SBarry Smith   pc->ops->view            = PCView_PFMG;
2465ebc551c0SBarry Smith   pc->ops->destroy         = PCDestroy_PFMG;
2466f91d8e95SBarry Smith   pc->ops->apply           = PCApply_PFMG;
24679e5bc791SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_PFMG;
246868326731SBarry Smith   pc->ops->setup           = PCSetUp_PFMG;
24692fa5cd67SKarl Rupp 
24709566063dSJacob Faibussowitsch   PetscCall(PetscCommGetComm(PetscObjectComm((PetscObject)pc), &ex->hcomm));
2471792fecdfSBarry Smith   PetscCallExternal(HYPRE_StructPFMGCreate, ex->hcomm, &ex->hsolver);
2472ebc551c0SBarry Smith   PetscFunctionReturn(0);
2473ebc551c0SBarry Smith }
2474d851a50bSGlenn Hammond 
2475d851a50bSGlenn Hammond /* we know we are working with a HYPRE_SStructMatrix */
2476d851a50bSGlenn Hammond typedef struct {
2477d851a50bSGlenn Hammond   MPI_Comm            hcomm; /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */
2478d851a50bSGlenn Hammond   HYPRE_SStructSolver ss_solver;
2479d851a50bSGlenn Hammond 
2480d851a50bSGlenn Hammond   /* keep copy of SYSPFMG options used so may view them */
24814ddd07fcSJed Brown   PetscInt its;
2482d851a50bSGlenn Hammond   double   tol;
24834ddd07fcSJed Brown   PetscInt relax_type;
24844ddd07fcSJed Brown   PetscInt num_pre_relax, num_post_relax;
2485d851a50bSGlenn Hammond } PC_SysPFMG;
2486d851a50bSGlenn Hammond 
24879371c9d4SSatish Balay PetscErrorCode PCDestroy_SysPFMG(PC pc) {
2488d851a50bSGlenn Hammond   PC_SysPFMG *ex = (PC_SysPFMG *)pc->data;
2489d851a50bSGlenn Hammond 
2490d851a50bSGlenn Hammond   PetscFunctionBegin;
2491792fecdfSBarry Smith   if (ex->ss_solver) PetscCallExternal(HYPRE_SStructSysPFMGDestroy, ex->ss_solver);
24929566063dSJacob Faibussowitsch   PetscCall(PetscCommRestoreComm(PetscObjectComm((PetscObject)pc), &ex->hcomm));
24939566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
2494d851a50bSGlenn Hammond   PetscFunctionReturn(0);
2495d851a50bSGlenn Hammond }
2496d851a50bSGlenn Hammond 
2497d851a50bSGlenn Hammond static const char *SysPFMGRelaxType[] = {"Weighted-Jacobi", "Red/Black-Gauss-Seidel"};
2498d851a50bSGlenn Hammond 
24999371c9d4SSatish Balay PetscErrorCode PCView_SysPFMG(PC pc, PetscViewer viewer) {
2500ace3abfcSBarry Smith   PetscBool   iascii;
2501d851a50bSGlenn Hammond   PC_SysPFMG *ex = (PC_SysPFMG *)pc->data;
2502d851a50bSGlenn Hammond 
2503d851a50bSGlenn Hammond   PetscFunctionBegin;
25049566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2505d851a50bSGlenn Hammond   if (iascii) {
25069566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  HYPRE SysPFMG preconditioning\n"));
250763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  max iterations %" PetscInt_FMT "\n", ex->its));
25089566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  tolerance %g\n", ex->tol));
25099566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  relax type %s\n", PFMGRelaxType[ex->relax_type]));
251063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  number pre-relax %" PetscInt_FMT " post-relax %" PetscInt_FMT "\n", ex->num_pre_relax, ex->num_post_relax));
2511d851a50bSGlenn Hammond   }
2512d851a50bSGlenn Hammond   PetscFunctionReturn(0);
2513d851a50bSGlenn Hammond }
2514d851a50bSGlenn Hammond 
25159371c9d4SSatish Balay PetscErrorCode PCSetFromOptions_SysPFMG(PC pc, PetscOptionItems *PetscOptionsObject) {
2516d851a50bSGlenn Hammond   PC_SysPFMG *ex  = (PC_SysPFMG *)pc->data;
2517ace3abfcSBarry Smith   PetscBool   flg = PETSC_FALSE;
2518d851a50bSGlenn Hammond 
2519d851a50bSGlenn Hammond   PetscFunctionBegin;
2520d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SysPFMG options");
25219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_syspfmg_print_statistics", "Print statistics", "HYPRE_SStructSysPFMGSetPrintLevel", flg, &flg, NULL));
252248a46eb9SPierre Jolivet   if (flg) PetscCallExternal(HYPRE_SStructSysPFMGSetPrintLevel, ex->ss_solver, 3);
25239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_syspfmg_its", "Number of iterations of SysPFMG to use as preconditioner", "HYPRE_SStructSysPFMGSetMaxIter", ex->its, &ex->its, NULL));
2524792fecdfSBarry Smith   PetscCallExternal(HYPRE_SStructSysPFMGSetMaxIter, ex->ss_solver, ex->its);
25259566063dSJacob Faibussowitsch   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));
2526792fecdfSBarry Smith   PetscCallExternal(HYPRE_SStructSysPFMGSetNumPreRelax, ex->ss_solver, ex->num_pre_relax);
25279566063dSJacob Faibussowitsch   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));
2528792fecdfSBarry Smith   PetscCallExternal(HYPRE_SStructSysPFMGSetNumPostRelax, ex->ss_solver, ex->num_post_relax);
2529d851a50bSGlenn Hammond 
25309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_syspfmg_tol", "Tolerance of SysPFMG", "HYPRE_SStructSysPFMGSetTol", ex->tol, &ex->tol, NULL));
2531792fecdfSBarry Smith   PetscCallExternal(HYPRE_SStructSysPFMGSetTol, ex->ss_solver, ex->tol);
2532dd39110bSPierre Jolivet   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));
2533792fecdfSBarry Smith   PetscCallExternal(HYPRE_SStructSysPFMGSetRelaxType, ex->ss_solver, ex->relax_type);
2534d0609cedSBarry Smith   PetscOptionsHeadEnd();
2535d851a50bSGlenn Hammond   PetscFunctionReturn(0);
2536d851a50bSGlenn Hammond }
2537d851a50bSGlenn Hammond 
25389371c9d4SSatish Balay PetscErrorCode PCApply_SysPFMG(PC pc, Vec x, Vec y) {
2539d851a50bSGlenn Hammond   PC_SysPFMG        *ex = (PC_SysPFMG *)pc->data;
2540d9ca1df4SBarry Smith   PetscScalar       *yy;
2541d9ca1df4SBarry Smith   const PetscScalar *xx;
25424ddd07fcSJed Brown   PetscInt           ilower[3], iupper[3];
25432cf14000SStefano Zampini   HYPRE_Int          hlower[3], hupper[3];
2544d851a50bSGlenn Hammond   Mat_HYPRESStruct  *mx       = (Mat_HYPRESStruct *)(pc->pmat->data);
25454ddd07fcSJed Brown   PetscInt           ordering = mx->dofs_order;
25464ddd07fcSJed Brown   PetscInt           nvars    = mx->nvars;
25474ddd07fcSJed Brown   PetscInt           part     = 0;
25484ddd07fcSJed Brown   PetscInt           size;
25494ddd07fcSJed Brown   PetscInt           i;
2550d851a50bSGlenn Hammond 
2551d851a50bSGlenn Hammond   PetscFunctionBegin;
25529566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(hypreCitation, &cite));
25539566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(mx->da, &ilower[0], &ilower[1], &ilower[2], &iupper[0], &iupper[1], &iupper[2]));
25542cf14000SStefano Zampini   /* when HYPRE_MIXEDINT is defined, sizeof(HYPRE_Int) == 32 */
2555d851a50bSGlenn Hammond   iupper[0] += ilower[0] - 1;
2556d851a50bSGlenn Hammond   iupper[1] += ilower[1] - 1;
2557d851a50bSGlenn Hammond   iupper[2] += ilower[2] - 1;
25582cf14000SStefano Zampini   hlower[0] = (HYPRE_Int)ilower[0];
25592cf14000SStefano Zampini   hlower[1] = (HYPRE_Int)ilower[1];
25602cf14000SStefano Zampini   hlower[2] = (HYPRE_Int)ilower[2];
25612cf14000SStefano Zampini   hupper[0] = (HYPRE_Int)iupper[0];
25622cf14000SStefano Zampini   hupper[1] = (HYPRE_Int)iupper[1];
25632cf14000SStefano Zampini   hupper[2] = (HYPRE_Int)iupper[2];
2564d851a50bSGlenn Hammond 
2565d851a50bSGlenn Hammond   size = 1;
25662fa5cd67SKarl Rupp   for (i = 0; i < 3; i++) size *= (iupper[i] - ilower[i] + 1);
25672fa5cd67SKarl Rupp 
2568d851a50bSGlenn Hammond   /* copy x values over to hypre for variable ordering */
2569d851a50bSGlenn Hammond   if (ordering) {
2570792fecdfSBarry Smith     PetscCallExternal(HYPRE_SStructVectorSetConstantValues, mx->ss_b, 0.0);
25719566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(x, &xx));
2572792fecdfSBarry Smith     for (i = 0; i < nvars; i++) PetscCallExternal(HYPRE_SStructVectorSetBoxValues, mx->ss_b, part, hlower, hupper, i, (HYPRE_Complex *)(xx + (size * i)));
25739566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(x, &xx));
2574792fecdfSBarry Smith     PetscCallExternal(HYPRE_SStructVectorAssemble, mx->ss_b);
2575792fecdfSBarry Smith     PetscCallExternal(HYPRE_SStructMatrixMatvec, 1.0, mx->ss_mat, mx->ss_b, 0.0, mx->ss_x);
2576792fecdfSBarry Smith     PetscCallExternal(HYPRE_SStructSysPFMGSolve, ex->ss_solver, mx->ss_mat, mx->ss_b, mx->ss_x);
2577d851a50bSGlenn Hammond 
2578d851a50bSGlenn Hammond     /* copy solution values back to PETSc */
25799566063dSJacob Faibussowitsch     PetscCall(VecGetArray(y, &yy));
2580792fecdfSBarry Smith     for (i = 0; i < nvars; i++) PetscCallExternal(HYPRE_SStructVectorGetBoxValues, mx->ss_x, part, hlower, hupper, i, (HYPRE_Complex *)(yy + (size * i)));
25819566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(y, &yy));
2582a65764d7SBarry Smith   } else { /* nodal ordering must be mapped to variable ordering for sys_pfmg */
2583d851a50bSGlenn Hammond     PetscScalar *z;
25844ddd07fcSJed Brown     PetscInt     j, k;
2585d851a50bSGlenn Hammond 
25869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nvars * size, &z));
2587792fecdfSBarry Smith     PetscCallExternal(HYPRE_SStructVectorSetConstantValues, mx->ss_b, 0.0);
25889566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(x, &xx));
2589d851a50bSGlenn Hammond 
2590d851a50bSGlenn Hammond     /* transform nodal to hypre's variable ordering for sys_pfmg */
2591d851a50bSGlenn Hammond     for (i = 0; i < size; i++) {
2592d851a50bSGlenn Hammond       k = i * nvars;
25932fa5cd67SKarl Rupp       for (j = 0; j < nvars; j++) z[j * size + i] = xx[k + j];
2594d851a50bSGlenn Hammond     }
2595792fecdfSBarry Smith     for (i = 0; i < nvars; i++) PetscCallExternal(HYPRE_SStructVectorSetBoxValues, mx->ss_b, part, hlower, hupper, i, (HYPRE_Complex *)(z + (size * i)));
25969566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(x, &xx));
2597792fecdfSBarry Smith     PetscCallExternal(HYPRE_SStructVectorAssemble, mx->ss_b);
2598792fecdfSBarry Smith     PetscCallExternal(HYPRE_SStructSysPFMGSolve, ex->ss_solver, mx->ss_mat, mx->ss_b, mx->ss_x);
2599d851a50bSGlenn Hammond 
2600d851a50bSGlenn Hammond     /* copy solution values back to PETSc */
26019566063dSJacob Faibussowitsch     PetscCall(VecGetArray(y, &yy));
2602792fecdfSBarry Smith     for (i = 0; i < nvars; i++) PetscCallExternal(HYPRE_SStructVectorGetBoxValues, mx->ss_x, part, hlower, hupper, i, (HYPRE_Complex *)(z + (size * i)));
2603d851a50bSGlenn Hammond     /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
2604d851a50bSGlenn Hammond     for (i = 0; i < size; i++) {
2605d851a50bSGlenn Hammond       k = i * nvars;
26062fa5cd67SKarl Rupp       for (j = 0; j < nvars; j++) yy[k + j] = z[j * size + i];
2607d851a50bSGlenn Hammond     }
26089566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(y, &yy));
26099566063dSJacob Faibussowitsch     PetscCall(PetscFree(z));
2610d851a50bSGlenn Hammond   }
2611d851a50bSGlenn Hammond   PetscFunctionReturn(0);
2612d851a50bSGlenn Hammond }
2613d851a50bSGlenn Hammond 
26149371c9d4SSatish Balay 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) {
2615d851a50bSGlenn Hammond   PC_SysPFMG *jac = (PC_SysPFMG *)pc->data;
26162cf14000SStefano Zampini   HYPRE_Int   oits;
2617d851a50bSGlenn Hammond 
2618d851a50bSGlenn Hammond   PetscFunctionBegin;
26199566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(hypreCitation, &cite));
2620792fecdfSBarry Smith   PetscCallExternal(HYPRE_SStructSysPFMGSetMaxIter, jac->ss_solver, its * jac->its);
2621792fecdfSBarry Smith   PetscCallExternal(HYPRE_SStructSysPFMGSetTol, jac->ss_solver, rtol);
26229566063dSJacob Faibussowitsch   PetscCall(PCApply_SysPFMG(pc, b, y));
2623792fecdfSBarry Smith   PetscCallExternal(HYPRE_SStructSysPFMGGetNumIterations, jac->ss_solver, &oits);
2624d851a50bSGlenn Hammond   *outits = oits;
2625d851a50bSGlenn Hammond   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
2626d851a50bSGlenn Hammond   else *reason = PCRICHARDSON_CONVERGED_RTOL;
2627792fecdfSBarry Smith   PetscCallExternal(HYPRE_SStructSysPFMGSetTol, jac->ss_solver, jac->tol);
2628792fecdfSBarry Smith   PetscCallExternal(HYPRE_SStructSysPFMGSetMaxIter, jac->ss_solver, jac->its);
2629d851a50bSGlenn Hammond   PetscFunctionReturn(0);
2630d851a50bSGlenn Hammond }
2631d851a50bSGlenn Hammond 
26329371c9d4SSatish Balay PetscErrorCode PCSetUp_SysPFMG(PC pc) {
2633d851a50bSGlenn Hammond   PC_SysPFMG       *ex = (PC_SysPFMG *)pc->data;
2634d851a50bSGlenn Hammond   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct *)(pc->pmat->data);
2635ace3abfcSBarry Smith   PetscBool         flg;
2636d851a50bSGlenn Hammond 
2637d851a50bSGlenn Hammond   PetscFunctionBegin;
26389566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATHYPRESSTRUCT, &flg));
263928b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Must use MATHYPRESSTRUCT with this preconditioner");
2640d851a50bSGlenn Hammond 
2641d851a50bSGlenn Hammond   /* create the hypre sstruct solver object and set its information */
2642792fecdfSBarry Smith   if (ex->ss_solver) PetscCallExternal(HYPRE_SStructSysPFMGDestroy, ex->ss_solver);
2643792fecdfSBarry Smith   PetscCallExternal(HYPRE_SStructSysPFMGCreate, ex->hcomm, &ex->ss_solver);
2644792fecdfSBarry Smith   PetscCallExternal(HYPRE_SStructSysPFMGSetZeroGuess, ex->ss_solver);
2645792fecdfSBarry Smith   PetscCallExternal(HYPRE_SStructSysPFMGSetup, ex->ss_solver, mx->ss_mat, mx->ss_b, mx->ss_x);
2646d851a50bSGlenn Hammond   PetscFunctionReturn(0);
2647d851a50bSGlenn Hammond }
2648d851a50bSGlenn Hammond 
2649d851a50bSGlenn Hammond /*MC
2650f1580f4eSBarry Smith      PCSYSPFMG - the hypre SysPFMG multigrid solver
2651d851a50bSGlenn Hammond 
2652d851a50bSGlenn Hammond    Level: advanced
2653d851a50bSGlenn Hammond 
2654f1580f4eSBarry Smith    Options Database Keys:
265567b8a455SSatish Balay + -pc_syspfmg_its <its> - number of iterations of SysPFMG to use as preconditioner
265667b8a455SSatish Balay . -pc_syspfmg_num_pre_relax <steps> - number of smoothing steps before coarse grid
265767b8a455SSatish Balay . -pc_syspfmg_num_post_relax <steps> - number of smoothing steps after coarse grid
265867b8a455SSatish Balay . -pc_syspfmg_tol <tol> - tolerance of SysPFMG
265967b8a455SSatish Balay - -pc_syspfmg_relax_type <Weighted-Jacobi,Red/Black-Gauss-Seidel> - relaxation type for the up and down cycles
2660d851a50bSGlenn Hammond 
266195452b02SPatrick Sanan    Notes:
2662f1580f4eSBarry Smith    See `PCPFMG` for hypre's PFMG that works for a scalar PDE and `PCSMG`
2663f1580f4eSBarry Smith 
2664f1580f4eSBarry Smith    See `PCHYPRE` for hypre's BoomerAMG algebraic multigrid solver
2665f1580f4eSBarry Smith 
266695452b02SPatrick Sanan    This is for CELL-centered descretizations
2667d851a50bSGlenn Hammond 
2668f1580f4eSBarry Smith    This must be used with the `MATHYPRESSTRUCT` matrix type.
2669d851a50bSGlenn Hammond 
2670f1580f4eSBarry Smith    This does not give access to all the functionality of hypres SysPFMG, it supports only one part, and one block per process defined by a PETSc `DMDA`.
2671f1580f4eSBarry Smith 
2672f1580f4eSBarry Smith .seealso: `PCMG`, `MATHYPRESSTRUCT`, `PCPFMG`, `PCHYPRE`, `PCGAMG`, `PCSMG`
2673d851a50bSGlenn Hammond M*/
2674d851a50bSGlenn Hammond 
26759371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_SysPFMG(PC pc) {
2676d851a50bSGlenn Hammond   PC_SysPFMG *ex;
2677d851a50bSGlenn Hammond 
2678d851a50bSGlenn Hammond   PetscFunctionBegin;
26799371c9d4SSatish Balay   PetscCall(PetscNew(&ex));
2680d851a50bSGlenn Hammond   pc->data = ex;
2681d851a50bSGlenn Hammond 
2682d851a50bSGlenn Hammond   ex->its            = 1;
2683d851a50bSGlenn Hammond   ex->tol            = 1.e-8;
2684d851a50bSGlenn Hammond   ex->relax_type     = 1;
2685d851a50bSGlenn Hammond   ex->num_pre_relax  = 1;
2686d851a50bSGlenn Hammond   ex->num_post_relax = 1;
2687d851a50bSGlenn Hammond 
2688d851a50bSGlenn Hammond   pc->ops->setfromoptions  = PCSetFromOptions_SysPFMG;
2689d851a50bSGlenn Hammond   pc->ops->view            = PCView_SysPFMG;
2690d851a50bSGlenn Hammond   pc->ops->destroy         = PCDestroy_SysPFMG;
2691d851a50bSGlenn Hammond   pc->ops->apply           = PCApply_SysPFMG;
2692d851a50bSGlenn Hammond   pc->ops->applyrichardson = PCApplyRichardson_SysPFMG;
2693d851a50bSGlenn Hammond   pc->ops->setup           = PCSetUp_SysPFMG;
26942fa5cd67SKarl Rupp 
26959566063dSJacob Faibussowitsch   PetscCall(PetscCommGetComm(PetscObjectComm((PetscObject)pc), &ex->hcomm));
2696792fecdfSBarry Smith   PetscCallExternal(HYPRE_SStructSysPFMGCreate, ex->hcomm, &ex->ss_solver);
2697d851a50bSGlenn Hammond   PetscFunctionReturn(0);
2698d851a50bSGlenn Hammond }
26991c188c59Sftrigaux 
2700f1580f4eSBarry Smith /* PC SMG */
27011c188c59Sftrigaux typedef struct {
27021c188c59Sftrigaux   MPI_Comm           hcomm; /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */
27031c188c59Sftrigaux   HYPRE_StructSolver hsolver;
27041c188c59Sftrigaux   PetscInt           its; /* keep copy of SMG options used so may view them */
27051c188c59Sftrigaux   double             tol;
27061c188c59Sftrigaux   PetscBool          print_statistics;
27071c188c59Sftrigaux   PetscInt           num_pre_relax, num_post_relax;
27081c188c59Sftrigaux } PC_SMG;
27091c188c59Sftrigaux 
27109371c9d4SSatish Balay PetscErrorCode PCDestroy_SMG(PC pc) {
27111c188c59Sftrigaux   PC_SMG *ex = (PC_SMG *)pc->data;
27121c188c59Sftrigaux 
27131c188c59Sftrigaux   PetscFunctionBegin;
27141c188c59Sftrigaux   if (ex->hsolver) PetscCallExternal(HYPRE_StructSMGDestroy, ex->hsolver);
27151c188c59Sftrigaux   PetscCall(PetscCommRestoreComm(PetscObjectComm((PetscObject)pc), &ex->hcomm));
27161c188c59Sftrigaux   PetscCall(PetscFree(pc->data));
27171c188c59Sftrigaux   PetscFunctionReturn(0);
27181c188c59Sftrigaux }
27191c188c59Sftrigaux 
27209371c9d4SSatish Balay PetscErrorCode PCView_SMG(PC pc, PetscViewer viewer) {
27211c188c59Sftrigaux   PetscBool iascii;
27221c188c59Sftrigaux   PC_SMG   *ex = (PC_SMG *)pc->data;
27231c188c59Sftrigaux 
27241c188c59Sftrigaux   PetscFunctionBegin;
27251c188c59Sftrigaux   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
27261c188c59Sftrigaux   if (iascii) {
27271c188c59Sftrigaux     PetscCall(PetscViewerASCIIPrintf(viewer, "  HYPRE SMG preconditioning\n"));
27281c188c59Sftrigaux     PetscCall(PetscViewerASCIIPrintf(viewer, "    max iterations %" PetscInt_FMT "\n", ex->its));
27291c188c59Sftrigaux     PetscCall(PetscViewerASCIIPrintf(viewer, "    tolerance %g\n", ex->tol));
27301c188c59Sftrigaux     PetscCall(PetscViewerASCIIPrintf(viewer, "    number pre-relax %" PetscInt_FMT " post-relax %" PetscInt_FMT "\n", ex->num_pre_relax, ex->num_post_relax));
27311c188c59Sftrigaux   }
27321c188c59Sftrigaux   PetscFunctionReturn(0);
27331c188c59Sftrigaux }
27341c188c59Sftrigaux 
27359371c9d4SSatish Balay PetscErrorCode PCSetFromOptions_SMG(PC pc, PetscOptionItems *PetscOptionsObject) {
27361c188c59Sftrigaux   PC_SMG *ex = (PC_SMG *)pc->data;
27371c188c59Sftrigaux 
27381c188c59Sftrigaux   PetscFunctionBegin;
27391c188c59Sftrigaux   PetscOptionsHeadBegin(PetscOptionsObject, "SMG options");
27401c188c59Sftrigaux 
27411c188c59Sftrigaux   PetscCall(PetscOptionsInt("-pc_smg_its", "Number of iterations of SMG to use as preconditioner", "HYPRE_StructSMGSetMaxIter", ex->its, &ex->its, NULL));
27421c188c59Sftrigaux   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));
27431c188c59Sftrigaux   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));
27441c188c59Sftrigaux   PetscCall(PetscOptionsReal("-pc_smg_tol", "Tolerance of SMG", "HYPRE_StructSMGSetTol", ex->tol, &ex->tol, NULL));
27451c188c59Sftrigaux 
27461c188c59Sftrigaux   PetscOptionsHeadEnd();
27471c188c59Sftrigaux   PetscFunctionReturn(0);
27481c188c59Sftrigaux }
27491c188c59Sftrigaux 
27509371c9d4SSatish Balay PetscErrorCode PCApply_SMG(PC pc, Vec x, Vec y) {
27511c188c59Sftrigaux   PC_SMG            *ex = (PC_SMG *)pc->data;
27521c188c59Sftrigaux   PetscScalar       *yy;
27531c188c59Sftrigaux   const PetscScalar *xx;
27541c188c59Sftrigaux   PetscInt           ilower[3], iupper[3];
27551c188c59Sftrigaux   HYPRE_Int          hlower[3], hupper[3];
27561c188c59Sftrigaux   Mat_HYPREStruct   *mx = (Mat_HYPREStruct *)(pc->pmat->data);
27571c188c59Sftrigaux 
27581c188c59Sftrigaux   PetscFunctionBegin;
27591c188c59Sftrigaux   PetscCall(PetscCitationsRegister(hypreCitation, &cite));
27601c188c59Sftrigaux   PetscCall(DMDAGetCorners(mx->da, &ilower[0], &ilower[1], &ilower[2], &iupper[0], &iupper[1], &iupper[2]));
27611c188c59Sftrigaux   /* when HYPRE_MIXEDINT is defined, sizeof(HYPRE_Int) == 32 */
27621c188c59Sftrigaux   iupper[0] += ilower[0] - 1;
27631c188c59Sftrigaux   iupper[1] += ilower[1] - 1;
27641c188c59Sftrigaux   iupper[2] += ilower[2] - 1;
27651c188c59Sftrigaux   hlower[0] = (HYPRE_Int)ilower[0];
27661c188c59Sftrigaux   hlower[1] = (HYPRE_Int)ilower[1];
27671c188c59Sftrigaux   hlower[2] = (HYPRE_Int)ilower[2];
27681c188c59Sftrigaux   hupper[0] = (HYPRE_Int)iupper[0];
27691c188c59Sftrigaux   hupper[1] = (HYPRE_Int)iupper[1];
27701c188c59Sftrigaux   hupper[2] = (HYPRE_Int)iupper[2];
27711c188c59Sftrigaux 
27721c188c59Sftrigaux   /* copy x values over to hypre */
27731c188c59Sftrigaux   PetscCallExternal(HYPRE_StructVectorSetConstantValues, mx->hb, 0.0);
27741c188c59Sftrigaux   PetscCall(VecGetArrayRead(x, &xx));
27751c188c59Sftrigaux   PetscCallExternal(HYPRE_StructVectorSetBoxValues, mx->hb, hlower, hupper, (HYPRE_Complex *)xx);
27761c188c59Sftrigaux   PetscCall(VecRestoreArrayRead(x, &xx));
27771c188c59Sftrigaux   PetscCallExternal(HYPRE_StructVectorAssemble, mx->hb);
27781c188c59Sftrigaux   PetscCallExternal(HYPRE_StructSMGSolve, ex->hsolver, mx->hmat, mx->hb, mx->hx);
27791c188c59Sftrigaux 
27801c188c59Sftrigaux   /* copy solution values back to PETSc */
27811c188c59Sftrigaux   PetscCall(VecGetArray(y, &yy));
27821c188c59Sftrigaux   PetscCallExternal(HYPRE_StructVectorGetBoxValues, mx->hx, hlower, hupper, (HYPRE_Complex *)yy);
27831c188c59Sftrigaux   PetscCall(VecRestoreArray(y, &yy));
27841c188c59Sftrigaux   PetscFunctionReturn(0);
27851c188c59Sftrigaux }
27861c188c59Sftrigaux 
27879371c9d4SSatish Balay 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) {
27881c188c59Sftrigaux   PC_SMG   *jac = (PC_SMG *)pc->data;
27891c188c59Sftrigaux   HYPRE_Int oits;
27901c188c59Sftrigaux 
27911c188c59Sftrigaux   PetscFunctionBegin;
27921c188c59Sftrigaux   PetscCall(PetscCitationsRegister(hypreCitation, &cite));
27931c188c59Sftrigaux   PetscCallExternal(HYPRE_StructSMGSetMaxIter, jac->hsolver, its * jac->its);
27941c188c59Sftrigaux   PetscCallExternal(HYPRE_StructSMGSetTol, jac->hsolver, rtol);
27951c188c59Sftrigaux 
27961c188c59Sftrigaux   PetscCall(PCApply_SMG(pc, b, y));
27971c188c59Sftrigaux   PetscCallExternal(HYPRE_StructSMGGetNumIterations, jac->hsolver, &oits);
27981c188c59Sftrigaux   *outits = oits;
27991c188c59Sftrigaux   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
28001c188c59Sftrigaux   else *reason = PCRICHARDSON_CONVERGED_RTOL;
28011c188c59Sftrigaux   PetscCallExternal(HYPRE_StructSMGSetTol, jac->hsolver, jac->tol);
28021c188c59Sftrigaux   PetscCallExternal(HYPRE_StructSMGSetMaxIter, jac->hsolver, jac->its);
28031c188c59Sftrigaux   PetscFunctionReturn(0);
28041c188c59Sftrigaux }
28051c188c59Sftrigaux 
28069371c9d4SSatish Balay PetscErrorCode PCSetUp_SMG(PC pc) {
28071c188c59Sftrigaux   PetscInt         i, dim;
28081c188c59Sftrigaux   PC_SMG          *ex = (PC_SMG *)pc->data;
28091c188c59Sftrigaux   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
28101c188c59Sftrigaux   PetscBool        flg;
28111c188c59Sftrigaux   DMBoundaryType   p[3];
28121c188c59Sftrigaux   PetscInt         M[3];
28131c188c59Sftrigaux 
28141c188c59Sftrigaux   PetscFunctionBegin;
28151c188c59Sftrigaux   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATHYPRESTRUCT, &flg));
28161c188c59Sftrigaux   PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Must use MATHYPRESTRUCT with this preconditioner");
28171c188c59Sftrigaux 
28181c188c59Sftrigaux   PetscCall(DMDAGetInfo(mx->da, &dim, &M[0], &M[1], &M[2], 0, 0, 0, 0, 0, &p[0], &p[1], &p[2], 0));
28191c188c59Sftrigaux   // Check if power of 2 in periodic directions
28201c188c59Sftrigaux   for (i = 0; i < dim; i++) {
28211c188c59Sftrigaux     if (((M[i] & (M[i] - 1)) != 0) && (p[i] == DM_BOUNDARY_PERIODIC)) {
28221c188c59Sftrigaux       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]);
28231c188c59Sftrigaux     }
28241c188c59Sftrigaux   }
28251c188c59Sftrigaux 
28261c188c59Sftrigaux   /* create the hypre solver object and set its information */
28271c188c59Sftrigaux   if (ex->hsolver) PetscCallExternal(HYPRE_StructSMGDestroy, (ex->hsolver));
28281c188c59Sftrigaux   PetscCallExternal(HYPRE_StructSMGCreate, ex->hcomm, &ex->hsolver);
28291c188c59Sftrigaux   // The hypre options must be set here and not in SetFromOptions because it is created here!
28301c188c59Sftrigaux   PetscCallExternal(HYPRE_StructSMGSetMaxIter, ex->hsolver, ex->its);
28311c188c59Sftrigaux   PetscCallExternal(HYPRE_StructSMGSetNumPreRelax, ex->hsolver, ex->num_pre_relax);
28321c188c59Sftrigaux   PetscCallExternal(HYPRE_StructSMGSetNumPostRelax, ex->hsolver, ex->num_post_relax);
28331c188c59Sftrigaux   PetscCallExternal(HYPRE_StructSMGSetTol, ex->hsolver, ex->tol);
28341c188c59Sftrigaux 
28351c188c59Sftrigaux   PetscCallExternal(HYPRE_StructSMGSetup, ex->hsolver, mx->hmat, mx->hb, mx->hx);
28361c188c59Sftrigaux   PetscCallExternal(HYPRE_StructSMGSetZeroGuess, ex->hsolver);
28371c188c59Sftrigaux   PetscFunctionReturn(0);
28381c188c59Sftrigaux }
28391c188c59Sftrigaux 
28401c188c59Sftrigaux /*MC
28415cb80ecdSBarry Smith      PCSMG - the hypre (structured grid) SMG multigrid solver
28421c188c59Sftrigaux 
28431c188c59Sftrigaux    Level: advanced
28441c188c59Sftrigaux 
2845f1580f4eSBarry Smith    Options Database Keys:
28465cb80ecdSBarry Smith + -pc_smg_its <its> - number of iterations of SMG to use as preconditioner
28475cb80ecdSBarry Smith . -pc_smg_num_pre_relax <steps> - number of smoothing steps before coarse grid
28485cb80ecdSBarry Smith . -pc_smg_num_post_relax <steps> - number of smoothing steps after coarse grid
28495cb80ecdSBarry Smith - -pc_smg_tol <tol> - tolerance of SMG
28501c188c59Sftrigaux 
28511c188c59Sftrigaux    Notes:
28521c188c59Sftrigaux    This is for CELL-centered descretizations
28531c188c59Sftrigaux 
28545cb80ecdSBarry Smith    This must be used with the `MATHYPRESTRUCT` `MatType`.
28551c188c59Sftrigaux 
2856f1580f4eSBarry Smith    This does not provide all the functionality of  hypre's SMG solver, it supports only one block per process defined by a PETSc `DMDA`.
2857f1580f4eSBarry Smith 
2858f1580f4eSBarry Smith    See `PCSYSPFMG`, `PCSMG`, `PCPFMG`, and `PCHYPRE` for access to hypre's other preconditioners
2859f1580f4eSBarry Smith 
2860f1580f4eSBarry Smith .seealso:  `PCMG`, `MATHYPRESTRUCT`, `PCPFMG`, `PCSYSPFMG`, `PCHYPRE`, `PCGAMG`
28611c188c59Sftrigaux M*/
28621c188c59Sftrigaux 
28639371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_SMG(PC pc) {
28641c188c59Sftrigaux   PC_SMG *ex;
28651c188c59Sftrigaux 
28661c188c59Sftrigaux   PetscFunctionBegin;
28679371c9d4SSatish Balay   PetscCall(PetscNew(&ex));
28681c188c59Sftrigaux   pc->data = ex;
28691c188c59Sftrigaux 
28701c188c59Sftrigaux   ex->its            = 1;
28711c188c59Sftrigaux   ex->tol            = 1.e-8;
28721c188c59Sftrigaux   ex->num_pre_relax  = 1;
28731c188c59Sftrigaux   ex->num_post_relax = 1;
28741c188c59Sftrigaux 
28751c188c59Sftrigaux   pc->ops->setfromoptions  = PCSetFromOptions_SMG;
28761c188c59Sftrigaux   pc->ops->view            = PCView_SMG;
28771c188c59Sftrigaux   pc->ops->destroy         = PCDestroy_SMG;
28781c188c59Sftrigaux   pc->ops->apply           = PCApply_SMG;
28791c188c59Sftrigaux   pc->ops->applyrichardson = PCApplyRichardson_SMG;
28801c188c59Sftrigaux   pc->ops->setup           = PCSetUp_SMG;
28811c188c59Sftrigaux 
28821c188c59Sftrigaux   PetscCall(PetscCommGetComm(PetscObjectComm((PetscObject)pc), &ex->hcomm));
28831c188c59Sftrigaux   PetscCallExternal(HYPRE_StructSMGCreate, ex->hcomm, &ex->hsolver);
28841c188c59Sftrigaux   PetscFunctionReturn(0);
28851c188c59Sftrigaux }
2886