Actual source code: hypre.c

  1: /*
  2:    Provides an interface to the LLNL package hypre
  3: */
  4: #include <petsc/private/petscimpl.h>
  5: #include <petscpc.h>

  7: /*@
  8:   PCHYPRESetDiscreteGradient - Set the discrete gradient matrix for `PCHYPRE` type of AMS or ADS

 10:   Collective

 12:   Input Parameters:
 13: + pc - the preconditioning context
 14: - G  - the discrete gradient

 16:   Level: intermediate

 18:   Notes:
 19:   `G` should have as many rows as the number of edges and as many columns as the number of vertices in the mesh

 21:   Each row of `G` has 2 nonzeros, with column indexes being the global indexes of edge's endpoints: matrix entries are +1 and -1 depending on edge orientation

 23:   Developer Note:
 24:   This automatically converts the matrix to `MATHYPRE` if it is not already of that type

 26: .seealso: [](ch_ksp), `PCHYPRE`, `PCHYPRESetDiscreteCurl()`
 27: @*/
 28: PetscErrorCode PCHYPRESetDiscreteGradient(PC pc, Mat G)
 29: {
 30:   PetscFunctionBegin;
 33:   PetscCheckSameComm(pc, 1, G, 2);
 34:   PetscTryMethod(pc, "PCHYPRESetDiscreteGradient_C", (PC, Mat), (pc, G));
 35:   PetscFunctionReturn(PETSC_SUCCESS);
 36: }

 38: /*@
 39:   PCHYPRESetDiscreteCurl - Set the discrete curl matrix for `PCHYPRE` type of ADS

 41:   Collective

 43:   Input Parameters:
 44: + pc - the preconditioning context
 45: - C  - the discrete curl

 47:   Level: intermediate

 49:   Notes:
 50:   `C` should have as many rows as the number of faces and as many columns as the number of edges in the mesh

 52:   Each row of `C` has as many nonzeros as the number of edges of a face, with column indexes being the global indexes of the corresponding edge.
 53:   Matrix entries are +1 and -1 depending on edge orientation with respect to the face orientation

 55:   Developer Notes:
 56:   This automatically converts the matrix to `MATHYPRE` if it is not already of that type

 58:   If this is only for  `PCHYPRE` type of ADS it should be called `PCHYPREADSSetDiscreteCurl()`

 60: .seealso: [](ch_ksp), `PCHYPRE`, `PCHYPRESetDiscreteGradient()`
 61: @*/
 62: PetscErrorCode PCHYPRESetDiscreteCurl(PC pc, Mat C)
 63: {
 64:   PetscFunctionBegin;
 67:   PetscCheckSameComm(pc, 1, C, 2);
 68:   PetscTryMethod(pc, "PCHYPRESetDiscreteCurl_C", (PC, Mat), (pc, C));
 69:   PetscFunctionReturn(PETSC_SUCCESS);
 70: }

 72: /*@
 73:   PCHYPRESetInterpolations - Set the interpolation matrices for `PCHYPRE` type of AMS or ADS

 75:   Collective

 77:   Input Parameters:
 78: + pc        - the preconditioning context
 79: . dim       - the dimension of the problem, only used in AMS
 80: . RT_PiFull - Raviart-Thomas interpolation matrix
 81: . RT_Pi     - x/y/z component of Raviart-Thomas interpolation matrix
 82: . ND_PiFull - Nedelec interpolation matrix
 83: - ND_Pi     - x/y/z component of Nedelec interpolation matrix

 85:   Level: intermediate

 87:   Notes:
 88:   For AMS, only Nedelec interpolation matrices are needed, the Raviart-Thomas interpolation matrices can be set to `NULL`.

 90:   For ADS, both type of interpolation matrices are needed.

 92:   Developer Note:
 93:   This automatically converts the matrix to `MATHYPRE` if it is not already of that type

 95: .seealso: [](ch_ksp), `PCHYPRE`
 96: @*/
 97: PetscErrorCode PCHYPRESetInterpolations(PC pc, PetscInt dim, Mat RT_PiFull, Mat RT_Pi[], Mat ND_PiFull, Mat ND_Pi[])
 98: {
 99:   PetscFunctionBegin;
101:   if (RT_PiFull) {
103:     PetscCheckSameComm(pc, 1, RT_PiFull, 3);
104:   }
105:   if (RT_Pi) {
106:     PetscAssertPointer(RT_Pi, 4);
107:     for (PetscInt i = 0; i < dim; ++i) {
108:       if (RT_Pi[i]) {
110:         PetscCheckSameComm(pc, 1, RT_Pi[i], 4);
111:       }
112:     }
113:   }
114:   if (ND_PiFull) {
116:     PetscCheckSameComm(pc, 1, ND_PiFull, 5);
117:   }
118:   if (ND_Pi) {
119:     PetscAssertPointer(ND_Pi, 6);
120:     for (PetscInt i = 0; i < dim; ++i) {
121:       if (ND_Pi[i]) {
123:         PetscCheckSameComm(pc, 1, ND_Pi[i], 6);
124:       }
125:     }
126:   }
127:   PetscTryMethod(pc, "PCHYPRESetInterpolations_C", (PC, PetscInt, Mat, Mat[], Mat, Mat[]), (pc, dim, RT_PiFull, RT_Pi, ND_PiFull, ND_Pi));
128:   PetscFunctionReturn(PETSC_SUCCESS);
129: }

131: /*@
132:   PCHYPRESetAlphaPoissonMatrix - Set the vector Poisson matrix for `PCHYPRE` of type AMS

134:   Collective

136:   Input Parameters:
137: + pc - the preconditioning context
138: - A  - the matrix

140:   Level: intermediate

142:   Note:
143:   `A` should be obtained by discretizing the vector valued Poisson problem with linear finite elements

145:   Developer Notes:
146:   This automatically converts the matrix to `MATHYPRE` if it is not already of that type

148:   If this is only for  `PCHYPRE` type of AMS it should be called `PCHYPREAMSSetAlphaPoissonMatrix()`

150: .seealso: [](ch_ksp), `PCHYPRE`, `PCHYPRESetDiscreteGradient()`, `PCHYPRESetDiscreteCurl()`, `PCHYPRESetBetaPoissonMatrix()`
151: @*/
152: PetscErrorCode PCHYPRESetAlphaPoissonMatrix(PC pc, Mat A)
153: {
154:   PetscFunctionBegin;
157:   PetscCheckSameComm(pc, 1, A, 2);
158:   PetscTryMethod(pc, "PCHYPRESetPoissonMatrix_C", (PC, Mat, PetscBool), (pc, A, PETSC_TRUE));
159:   PetscFunctionReturn(PETSC_SUCCESS);
160: }

162: /*@
163:   PCHYPRESetBetaPoissonMatrix - Set the Poisson matrix for `PCHYPRE` of type AMS

165:   Collective

167:   Input Parameters:
168: + pc - the preconditioning context
169: - A  - the matrix, or `NULL` to turn it off

171:   Level: intermediate

173:   Note:
174:   `A` should be obtained by discretizing the Poisson problem with linear finite elements.

176:   Developer Notes:
177:   This automatically converts the matrix to `MATHYPRE` if it is not already of that type

179:   If this is only for  `PCHYPRE` type of AMS it should be called `PCHYPREAMSPCHYPRESetBetaPoissonMatrix()`

181: .seealso: [](ch_ksp), `PCHYPRE`, `PCHYPRESetDiscreteGradient()`, `PCHYPRESetDiscreteCurl()`, `PCHYPRESetAlphaPoissonMatrix()`
182: @*/
183: PetscErrorCode PCHYPRESetBetaPoissonMatrix(PC pc, Mat A)
184: {
185:   PetscFunctionBegin;
187:   if (A) {
189:     PetscCheckSameComm(pc, 1, A, 2);
190:   }
191:   PetscTryMethod(pc, "PCHYPRESetPoissonMatrix_C", (PC, Mat, PetscBool), (pc, A, PETSC_FALSE));
192:   PetscFunctionReturn(PETSC_SUCCESS);
193: }

195: /*@
196:   PCHYPRESetEdgeConstantVectors - Set the representation of the constant vector fields in the edge element basis for `PCHYPRE` of type AMS

198:   Collective

200:   Input Parameters:
201: + pc  - the preconditioning context
202: . ozz - vector representing (1,0,0) (or (1,0) in 2D)
203: . zoz - vector representing (0,1,0) (or (0,1) in 2D)
204: - zzo - vector representing (0,0,1) (use NULL in 2D)

206:   Level: intermediate

208:   Developer Note:
209:   If this is only for  `PCHYPRE` type of AMS it should be called `PCHYPREAMSSetEdgeConstantVectors()`

211: .seealso: [](ch_ksp), `PCHYPRE`, `PCHYPRESetDiscreteGradient()`, `PCHYPRESetDiscreteCurl()`, `PCHYPRESetAlphaPoissonMatrix()`
212: @*/
213: PetscErrorCode PCHYPRESetEdgeConstantVectors(PC pc, Vec ozz, Vec zoz, Vec zzo)
214: {
215:   PetscFunctionBegin;
220:   PetscCheckSameComm(pc, 1, ozz, 2);
221:   PetscCheckSameComm(pc, 1, zoz, 3);
222:   if (zzo) PetscCheckSameComm(pc, 1, zzo, 4);
223:   PetscTryMethod(pc, "PCHYPRESetEdgeConstantVectors_C", (PC, Vec, Vec, Vec), (pc, ozz, zoz, zzo));
224:   PetscFunctionReturn(PETSC_SUCCESS);
225: }

227: /*@
228:   PCHYPREAMSSetInteriorNodes - Set the list of interior nodes to a zero-conductivity region for `PCHYPRE` of type AMS

230:   Collective

232:   Input Parameters:
233: + pc       - the preconditioning context
234: - interior - vector. node is interior if its entry in the array is 1.0.

236:   Level: intermediate

238:   Note:
239:   This calls `HYPRE_AMSSetInteriorNodes()`

241: .seealso: [](ch_ksp), `PCHYPRE`, `PCHYPRESetDiscreteGradient()`, `PCHYPRESetDiscreteCurl()`, `PCHYPRESetAlphaPoissonMatrix()`
242: @*/
243: PetscErrorCode PCHYPREAMSSetInteriorNodes(PC pc, Vec interior)
244: {
245:   PetscFunctionBegin;
248:   PetscCheckSameComm(pc, 1, interior, 2);
249:   PetscTryMethod(pc, "PCHYPREAMSSetInteriorNodes_C", (PC, Vec), (pc, interior));
250:   PetscFunctionReturn(PETSC_SUCCESS);
251: }

253: /*@
254:   PCHYPRESetType - Sets which hypre preconditioner you wish to use

256:   Input Parameters:
257: + pc   - the preconditioner context
258: - name - either euclid, ilu, pilut, parasails, boomeramg, ams, or ads

260:   Options Database Key:
261: . pc_hypre_type - One of euclid, ilu, pilut, parasails, boomeramg, ams, or ads

263:   Level: intermediate

265: .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCHYPRE`
266: @*/
267: PetscErrorCode PCHYPRESetType(PC pc, const char name[])
268: {
269:   PetscFunctionBegin;
271:   PetscAssertPointer(name, 2);
272:   PetscTryMethod(pc, "PCHYPRESetType_C", (PC, const char[]), (pc, name));
273:   PetscFunctionReturn(PETSC_SUCCESS);
274: }

276: /*@C
277:   PCHYPREGetCFMarkers - Gets CF marker arrays for all levels (except the finest level)

279:   Logically Collective

281:   Input Parameter:
282: . pc - the preconditioner context

284:   Output Parameters:
285: + n_per_level - the number of nodes per level (size of `num_levels`)
286: - CFMarkers   - the Coarse/Fine Boolean arrays (size of `num_levels` - 1)

288:   Level: advanced

290:   Note:
291:   Caller is responsible for memory management of `n_per_level` and `CFMarkers` pointers. That is they should free them with `PetscFree()` when no longer needed.

293: .seealso: [](ch_ksp), `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()`
294: @*/
295: PetscErrorCode PCHYPREGetCFMarkers(PC pc, PetscInt *n_per_level[], PetscBT *CFMarkers[])
296: {
297:   PetscFunctionBegin;
299:   PetscAssertPointer(n_per_level, 2);
300:   PetscAssertPointer(CFMarkers, 3);
301:   PetscUseMethod(pc, "PCHYPREGetCFMarkers_C", (PC, PetscInt *[], PetscBT *[]), (pc, n_per_level, CFMarkers));
302:   PetscFunctionReturn(PETSC_SUCCESS);
303: }

305: /*@
306:   PCHYPREGetType - Gets which hypre preconditioner you are using

308:   Input Parameter:
309: . pc - the preconditioner context

311:   Output Parameter:
312: . name - either euclid, ilu, pilut, parasails, boomeramg, ams, or ads

314:   Level: intermediate

316: .seealso: [](ch_ksp), `PCCreate()`, `PCHYPRESetType()`, `PCType`, `PC`, `PCHYPRE`
317: @*/
318: PetscErrorCode PCHYPREGetType(PC pc, const char *name[])
319: {
320:   PetscFunctionBegin;
322:   PetscAssertPointer(name, 2);
323:   PetscTryMethod(pc, "PCHYPREGetType_C", (PC, const char *[]), (pc, name));
324:   PetscFunctionReturn(PETSC_SUCCESS);
325: }

327: /*@
328:   PCMGGalerkinSetMatProductAlgorithm - Set type of sparse matrix-matrix product for hypre's BoomerAMG to use on GPUs

330:   Logically Collective

332:   Input Parameters:
333: + pc   - the hypre context
334: - name - one of 'cusparse', 'hypre'

336:   Options Database Key:
337: . -pc_mg_galerkin_mat_product_algorithm (cusparse|hypre) - Type of sparse matrix-matrix product to use in hypre

339:   Level: intermediate

341:   Developer Note:
342:   How the name starts with `PCMG`, should it not be `PCHYPREBoomerAMG`?

344: .seealso: [](ch_ksp), `PCHYPRE`, `PCMGGalerkinGetMatProductAlgorithm()`
345: @*/
346: PetscErrorCode PCMGGalerkinSetMatProductAlgorithm(PC pc, const char name[])
347: {
348:   PetscFunctionBegin;
350:   PetscTryMethod(pc, "PCMGGalerkinSetMatProductAlgorithm_C", (PC, const char[]), (pc, name));
351:   PetscFunctionReturn(PETSC_SUCCESS);
352: }

354: /*@
355:   PCMGGalerkinGetMatProductAlgorithm - Get type of sparse matrix-matrix product for hypre's BoomerAMG to use on GPUs

357:   Not Collective

359:   Input Parameter:
360: . pc - the multigrid context

362:   Output Parameter:
363: . name - one of 'cusparse', 'hypre'

365:   Level: intermediate

367: .seealso: [](ch_ksp), `PCHYPRE`, `PCMGGalerkinSetMatProductAlgorithm()`
368: @*/
369: PetscErrorCode PCMGGalerkinGetMatProductAlgorithm(PC pc, const char *name[])
370: {
371:   PetscFunctionBegin;
373:   PetscTryMethod(pc, "PCMGGalerkinGetMatProductAlgorithm_C", (PC, const char *[]), (pc, name));
374:   PetscFunctionReturn(PETSC_SUCCESS);
375: }