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