Actual source code: hyppilut.c
1: /*$Id: bvec2.c,v 1.202 2001/09/12 03:26:24 bsmith Exp $*/
2: /*
4: */
6: #include src/sles/pc/pcimpl.h
7: EXTERN_C_BEGIN
8: #include "HYPRE.h"
9: #include "IJ_mv.h"
10: #include "HYPRE_parcsr_ls.h"
11: EXTERN_C_END
13: extern int MatHYPRE_IJMatrixCreate(Mat,HYPRE_IJMatrix*);
14: extern int MatHYPRE_IJMatrixCopy(Mat,HYPRE_IJMatrix);
15: extern int VecHYPRE_IJVectorCreate(Vec,HYPRE_IJVector*);
17: /*
18: Private context (data structure) for the preconditioner.
19: */
20: typedef struct {
21: HYPRE_Solver hsolver;
22: HYPRE_IJMatrix ij;
23: HYPRE_IJVector b,x;
25: int (*destroy)(HYPRE_Solver);
26: int (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
27: int (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
29: /* options for pilut and BoomerAMG*/
30: int maxiter;
31: double tol;
33: /* options for pilut */
34: int factorrowsize;
36: /* options for parasails */
37: int nlevels;
38: double threshhold;
39: double filter;
40: int sym;
41: double loadbal;
42: int logging;
43: int ruse;
44: int symt;
46: /* options for euclid */
47: PetscTruth bjilu;
48: int levels;
50: /* options for euclid and BoomerAMG */
51: PetscTruth printstatistics;
53: /* options for BoomerAMG */
54: int maxlevels;
55: double strongthreshold;
56: double maxrowsum;
57: int *gridsweeps;
58: int coarsentype;
59: int measuretype;
60: int *relaxtype;
61: int **gridrelaxpoints;
62: } PC_HYPRE;
65: static int PCSetUp_HYPRE(PC pc)
66: {
67: PC_HYPRE *jac = (PC_HYPRE*)pc->data;
68: int ierr;
69: HYPRE_ParCSRMatrix hmat;
70: HYPRE_ParVector bv,xv;
73: if (!jac->ij) { /* create the matrix the first time through */
74: MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);
75: }
76: if (!jac->b) {
77: VecHYPRE_IJVectorCreate(pc->vec,&jac->b);
78: VecHYPRE_IJVectorCreate(pc->vec,&jac->x);
79: }
80: MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);
81: HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);
82: HYPRE_IJVectorGetObject(jac->b,(void**)&bv);
83: HYPRE_IJVectorGetObject(jac->x,(void**)&xv);
84: (*jac->setup)(jac->hsolver,hmat,bv,xv);
85: return(0);
86: }
88: /*
89: Replaces the address where the HYPRE vector points to its data with the address of
90: PETSc's data. Saves the old address so it can be reset when we are finished with it.
91: Allows use to get the data into a HYPRE vector without the cost of memcopies
92: */
93: #define HYPREReplacePointer(b,newvalue,savedvalue) {
94: hypre_ParVector *par_vector = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));
95: hypre_Vector *local_vector = hypre_ParVectorLocalVector(par_vector);
96: savedvalue = local_vector->data;
97: local_vector->data = newvalue;}
99: static int PCApply_HYPRE(PC pc,Vec b,Vec x)
100: {
101: PC_HYPRE *jac = (PC_HYPRE*)pc->data;
102: int ierr;
103: HYPRE_ParCSRMatrix hmat;
104: PetscScalar *bv,*xv;
105: HYPRE_ParVector jbv,jxv;
106: PetscScalar *sbv,*sxv;
109: VecGetArray(b,&bv);
110: VecGetArray(x,&xv);
111: HYPREReplacePointer(jac->b,bv,sbv);
112: HYPREReplacePointer(jac->x,xv,sxv);
114: HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);
115: HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);
116: HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);
117: (*jac->solve)(jac->hsolver,hmat,jbv,jxv);
119: HYPREReplacePointer(jac->b,sbv,bv);
120: HYPREReplacePointer(jac->x,sxv,xv);
121: VecRestoreArray(x,&xv);
122: VecRestoreArray(b,&bv);
124: return(0);
125: }
127: static int PCDestroy_HYPRE(PC pc)
128: {
129: PC_HYPRE *jac = (PC_HYPRE*)pc->data;
130: int ierr;
133: HYPRE_IJMatrixDestroy(jac->ij);
134: HYPRE_IJVectorDestroy(jac->b);
135: HYPRE_IJVectorDestroy(jac->x);
136: (*jac->destroy)(jac->hsolver);
137: PetscFree(jac);
138: return(0);
139: }
141: /* --------------------------------------------------------------------------------------------*/
142: static int PCSetFromOptions_HYPRE_Pilut(PC pc)
143: {
144: PC_HYPRE *jac = (PC_HYPRE*)pc->data;
145: int ierr;
146: PetscTruth flag;
149: PetscOptionsHead("HYPRE Pilut Options");
150: PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);
151: if (flag) {
152: HYPRE_ParCSRPilutSetMaxIter(jac->hsolver,jac->maxiter);
153: }
154: PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);
155: if (flag) {
156: HYPRE_ParCSRPilutSetDropTolerance(jac->hsolver,jac->tol);
157: }
158: PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);
159: if (flag) {
160: HYPRE_ParCSRPilutSetFactorRowSize(jac->hsolver,jac->factorrowsize);
161: }
162: PetscOptionsTail();
163: return(0);
164: }
166: static int PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
167: {
168: PC_HYPRE *jac = (PC_HYPRE*)pc->data;
169: int ierr;
170: PetscTruth isascii;
173: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
174: if (isascii) {
175: PetscViewerASCIIPrintf(viewer," HYPRE Pilut preconditioningn");
176: if (jac->maxiter != PETSC_DEFAULT) {
177: PetscViewerASCIIPrintf(viewer," HYPRE Pilut: maximum number of iterations %dn",jac->maxiter);
178: } else {
179: PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default maximum number of iterations n");
180: }
181: if (jac->tol != PETSC_DEFAULT) {
182: PetscViewerASCIIPrintf(viewer," HYPRE Pilut: drop tolerance %gn",jac->tol);
183: } else {
184: PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default drop tolerance n");
185: }
186: if (jac->factorrowsize != PETSC_DEFAULT) {
187: PetscViewerASCIIPrintf(viewer," HYPRE Pilut: factor row size %dn",jac->factorrowsize);
188: } else {
189: PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default factor row size n");
190: }
191: }
192: return(0);
193: }
195: /* --------------------------------------------------------------------------------------------*/
196: static int PCSetFromOptions_HYPRE_Euclid(PC pc)
197: {
198: PC_HYPRE *jac = (PC_HYPRE*)pc->data;
199: int ierr;
200: PetscTruth flag;
201: char *args[2];
204: jac->bjilu = PETSC_FALSE;
205: jac->levels = 1;
207: PetscOptionsHead("HYPRE Euclid Options");
208: PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);
209: if (flag) {
210: char levels[16];
211: if (jac->levels < 0) SETERRQ1(1,"Number of levels %d must be nonegative",jac->levels);
212: sprintf(levels,"%d",jac->levels);
213: args[0] = "-level"; args[1] = levels;
214: HYPRE_EuclidSetParams(jac->hsolver,2,args);
215: }
216: PetscOptionsLogical("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);
217: if (jac->bjilu) {
218: args[0] = "-bj"; args[1] = "1";
219: HYPRE_EuclidSetParams(jac->hsolver,2,args);
220: }
221:
222: PetscOptionsLogical("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);
223: if (jac->printstatistics) {
224: args[0] = "-eu_stats"; args[1] = "1";
225: HYPRE_EuclidSetParams(jac->hsolver,2,args);
226: args[0] = "-eu_mem"; args[1] = "1";
227: HYPRE_EuclidSetParams(jac->hsolver,2,args);
228: }
229: PetscOptionsTail();
230: return(0);
231: }
233: static int PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
234: {
235: PC_HYPRE *jac = (PC_HYPRE*)pc->data;
236: int ierr;
237: PetscTruth isascii;
240: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
241: if (isascii) {
242: PetscViewerASCIIPrintf(viewer," HYPRE Euclid preconditioningn");
243: PetscViewerASCIIPrintf(viewer," HYPRE Euclid: number of levels %dn",jac->levels);
244: if (jac->bjilu) {
245: PetscViewerASCIIPrintf(viewer," HYPRE Euclid: Using block Jacobi ILU instead of parallel ILUn");
246: }
247: }
248: return(0);
249: }
251: /* --------------------------------------------------------------------------------------------*/
253: static char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout"};
254: static char *HYPREBoomerAMGMeasureType[] = {"local","global"};
255: static char *HYPREBoomerAMGRelaxType[] = {"Jacobi","sequential-Gauss-Seidel","","Gauss-Seidel/Jacobi","","","symmetric-Gauss-Seidel/Jacobi",
256: "","","Gaussian-elimination"};
257: static int PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
258: {
259: PC_HYPRE *jac = (PC_HYPRE*)pc->data;
260: int ierr,n = 4,i;
261: PetscTruth flg;
262: char result[32];
265: jac->maxlevels = 25;
266: jac->maxiter = 20;
267: jac->tol = 1.e-7;
268: jac->strongthreshold = .25;
269: jac->maxrowsum = .9;
270: jac->coarsentype = 6;
271: jac->measuretype = 0;
272:
274: /* this is terrible; HYPRE frees this array so we have to malloc it */
275: jac->gridsweeps = (int*)malloc(4*sizeof(int));
276: jac->gridsweeps[0] = jac->gridsweeps[1] = jac->gridsweeps[2] = 2;
277: jac->gridsweeps[3] = 1;
279: jac->relaxtype = (int*)malloc(4*sizeof(int));
280: jac->relaxtype[0] = jac->relaxtype[1] = jac->relaxtype[2] = 3;
281: jac->relaxtype[3] = 9;
283: PetscOptionsHead("HYPRE BoomerAMG Options");
284: PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);
285: if (flg) {
286: if (jac->maxlevels < 2) SETERRQ1(1,"Number of levels %d must be at least one",jac->maxlevels);
287: HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);
288: }
289: PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used","None",jac->maxiter,&jac->maxiter,&flg);
290: if (flg) {
291: if (jac->maxiter < 1) SETERRQ1(1,"Number of iterations %d must be at least one",jac->maxiter);
292: HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);
293: }
294: PetscOptionsScalar("-pc_hypre_boomeramg_tol","Convergence tolerance","None",jac->tol,&jac->tol,&flg);
295: if (flg) {
296: if (jac->tol < 0.0) SETERRQ1(1,"Tolerance %g must be great than or equal zero",jac->tol);
297: HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);
298: }
299: PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);
300: if (flg) {
301: if (jac->strongthreshold < 0.0) SETERRQ1(1,"Strong threshold %g must be great than or equal zero",jac->strongthreshold);
302: HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);
303: }
304: PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);
305: if (flg) {
306: if (jac->maxrowsum < 0.0) SETERRQ1(1,"Maximum row sum %g must be greater than zero",jac->maxrowsum);
307: if (jac->maxrowsum > 1.0) SETERRQ1(1,"Maximum row sum %g must be less than or equal one",jac->maxrowsum);
308: HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);
309: }
310:
311: n = 4;
312: PetscOptionsIntArray("-pc_hypre_boomeramg_grid_sweeps","Grid sweeps for fine,down,up,coarse","None",jac->gridsweeps,&n,&flg);
313: if (flg) {
314: if (n == 1) {
315: jac->gridsweeps[1] = jac->gridsweeps[2] = jac->gridsweeps[3] = jac->gridsweeps[0];
316: n = 4;
317: }
318: if (n != 4) SETERRQ1(1,"You must provide either 1 or 4 values seperated by commas, you provided %d",n);
319: HYPRE_BoomerAMGSetNumGridSweeps(jac->hsolver,jac->gridsweeps);
320: CHKMEMQ;
321: }
322: /*
323: Suggested by QUANDALLE Philippe <Philippe.QUANDALLE@ifp.fr>
325: gridrelaxpoints[i][j] are for i=0,1,2,3 (fine,up,down,coarse) and j=sweep number
326: 0 indicates smooth all points
327: 1 indicates smooth coarse points
328: -1 indicates smooth fine points
330: Here when j=1 it first smooths all the coarse points, then all the fine points.
331: */
332: jac->gridrelaxpoints = (int**)malloc(4*sizeof(int*));
333: if(jac->gridsweeps[0]>0) jac->gridrelaxpoints[0] = (int*)malloc(jac->gridsweeps[0]*sizeof(int));
334: if(jac->gridsweeps[1]>0) jac->gridrelaxpoints[1] = (int*)malloc(jac->gridsweeps[1]*sizeof(int));
335: if(jac->gridsweeps[2]>0) jac->gridrelaxpoints[2] = (int*)malloc(jac->gridsweeps[2]*sizeof(int));
336: if(jac->gridsweeps[3]>0) jac->gridrelaxpoints[3] = (int*)malloc(jac->gridsweeps[3]*sizeof(int));
338: PetscOptionsLogical("-pc_hypre_boomeramg_sweep_all","Sweep all points","None",PETSC_FALSE,&flg,0);
339: if(jac->gridsweeps[0] == 1) jac->gridrelaxpoints[0][0] = 0;
340: else if(jac->gridsweeps[0] == 2) {
341: if (flg) {
342: jac->gridrelaxpoints[0][0] = 0; jac->gridrelaxpoints[0][1] = 0;
343: } else {
344: jac->gridrelaxpoints[0][0] = 1; jac->gridrelaxpoints[0][1] = -1;
345: }
346: } else if (jac->gridsweeps[0] > 2) {
347: SETERRQ(1,"Grid sweeps can only be 0, 1, or 2");
348: }
349:
350: if(jac->gridsweeps[1] == 1) jac->gridrelaxpoints[1][0] = 0;
351: else if(jac->gridsweeps[1] == 2) {
352: if (flg) {
353: jac->gridrelaxpoints[1][0] = 0; jac->gridrelaxpoints[1][1] = 0;
354: } else {
355: jac->gridrelaxpoints[1][0] = -1; jac->gridrelaxpoints[1][1] = 1;
356: }
357: } else if (jac->gridsweeps[1] > 2) {
358: SETERRQ(1,"Grid sweeps can only be 0, 1, or 2");
359: }
361: if(jac->gridsweeps[2] == 1) jac->gridrelaxpoints[2][0] = 0;
362: else if(jac->gridsweeps[2] == 2) {
363: if (flg) {
364: jac->gridrelaxpoints[2][0] = 0; jac->gridrelaxpoints[2][1] = 0;
365: } else {
366: jac->gridrelaxpoints[2][0] = 1; jac->gridrelaxpoints[2][1] = -1;
367: }
368: } else if (jac->gridsweeps[2] > 2) {
369: SETERRQ(1,"Grid sweeps can only be 0, 1, or 2");
370: }
372: for (i=0; i<jac->gridsweeps[3]; i++) {
373: jac->gridrelaxpoints[3][i] = 0;
374: }
375: HYPRE_BoomerAMGSetGridRelaxPoints(jac->hsolver,jac->gridrelaxpoints);
378: PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],result,16,&flg);
379: if (flg) {
380: int type = -1;
381: for (i=0; i<2; i++) {
382: PetscStrcmp(result,HYPREBoomerAMGMeasureType[i],&flg);
383: if (flg) {
384: type = i;
385: break;
386: }
387: }
388: if (type == -1) SETERRQ1(1,"Unknown measure type %s",result);
389: jac->measuretype = type;
390: HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);
391: }
392: PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,7,HYPREBoomerAMGCoarsenType[6],result,25,&flg);
393: if (flg) {
394: int type = -1;
395: for (i=0; i<7; i++) {
396: PetscStrcmp(result,HYPREBoomerAMGCoarsenType[i],&flg);
397: if (flg) {
398: type = i;
399: break;
400: }
401: }
402: if (type == -1) SETERRQ1(1,"Unknown coarsen type %s",result);
403: jac->coarsentype = type;
404: HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);
405: }
406: PetscOptionsEList("-pc_hypre_boomeramg_relax_type","Relax type","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],result,32,&flg);
407: if (flg) {
408: int type = -1;
409: for (i=0; i<10; i++) {
410: PetscStrcmp(result,HYPREBoomerAMGRelaxType[i],&flg);
411: if (flg) {
412: type = i;
413: break;
414: }
415: }
416: if (type == -1) SETERRQ1(1,"Unknown relax type %s",result);
417: jac->relaxtype[0] = jac->relaxtype[1] = jac->relaxtype[2] = type;
418: }
419: PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],result,32,&flg);
420: if (flg) {
421: int type = -1;
422: for (i=0; i<10; i++) {
423: PetscStrcmp(result,HYPREBoomerAMGRelaxType[i],&flg);
424: if (flg) {
425: type = i;
426: break;
427: }
428: }
429: if (type == -1) SETERRQ1(1,"Unknown relax type %s",result);
430: jac->relaxtype[3] = type;
431: }
432: HYPRE_BoomerAMGSetGridRelaxType(jac->hsolver,jac->relaxtype);
433:
434: PetscOptionsLogical("-pc_hypre_boomeramg_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);
435: if (jac->printstatistics) {
436: HYPRE_BoomerAMGSetIOutDat(jac->hsolver,3);
437: }
438: PetscOptionsTail();
439: return(0);
440: }
442: static int PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
443: {
444: PC_HYPRE *jac = (PC_HYPRE*)pc->data;
445: int ierr;
446: PetscTruth isascii;
449: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
450: if (isascii) {
451: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG preconditioningn");
452: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum number of levels %dn",jac->maxlevels);
453: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum number of iterations %dn",jac->maxiter);
454: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Convergence tolerance %gn",jac->tol);
455: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Threshold for strong coupling %gn",jac->strongthreshold);
456: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum row sums %gn",jac->maxrowsum);
457: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps on fine grid %dn",jac->gridsweeps[0]);
458: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps down %dn",jac->gridsweeps[1]);
459: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps up %dn",jac->gridsweeps[2]);
460: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps on coarse %dn",jac->gridsweeps[3]);
462: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax on fine grid %sn",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);
463: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax down %sn",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);
464: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax up %sn",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);
465: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax on coarse %sn",HYPREBoomerAMGRelaxType[jac->relaxtype[3]]);
467: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Measure type %sn",HYPREBoomerAMGMeasureType[jac->measuretype]);
468: PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Coarsen type %sn",HYPREBoomerAMGCoarsenType[jac->coarsentype]);
469: }
470: return(0);
471: }
473: /* --------------------------------------------------------------------------------------------*/
474: static int PCSetFromOptions_HYPRE_ParaSails(PC pc)
475: {
476: PC_HYPRE *jac = (PC_HYPRE*)pc->data;
477: int ierr;
478: PetscTruth flag;
479: char *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"},buff[32];
482: jac->nlevels = 1;
483: jac->threshhold = .1;
484: jac->filter = .1;
485: jac->loadbal = 0;
486: if (PetscLogPrintInfo) {
487: jac->logging = (int) PETSC_TRUE;
488: } else {
489: jac->logging = (int) PETSC_FALSE;
490: }
491: jac->ruse = (int) PETSC_TRUE;
492: jac->symt = 0;
494: PetscOptionsHead("HYPRE ParaSails Options");
495: PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);
496: PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,0);
497: HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);
499: PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,0);
500: HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);
502: PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,0);
503: HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);
505: PetscOptionsLogical("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,0);
506: HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);
508: PetscOptionsLogical("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,0);
509: HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);
511: PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],buff,32,&flag);
512: while (flag) {
513: PetscStrcmp(symtlist[0],buff,&flag);
514: if (flag) {
515: jac->symt = 0;
516: break;
517: }
518: PetscStrcmp(symtlist[1],buff,&flag);
519: if (flag) {
520: jac->symt = 1;
521: break;
522: }
523: PetscStrcmp(symtlist[2],buff,&flag);
524: if (flag) {
525: jac->symt = 2;
526: break;
527: }
528: SETERRQ1(1,"Unknown HYPRE ParaSails Sym option %s",buff);
529: }
530: HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);
532: PetscOptionsTail();
533: return(0);
534: }
536: static int PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
537: {
538: PC_HYPRE *jac = (PC_HYPRE*)pc->data;
539: int ierr;
540: PetscTruth isascii;
541: char *symt;
544: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
545: if (isascii) {
546: PetscViewerASCIIPrintf(viewer," HYPRE ParaSails preconditioningn");
547: PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: nlevels %dn",jac->nlevels);
548: PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: threshold %gn",jac->threshhold);
549: PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: filter %gn",jac->filter);
550: PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: load balance %gn",jac->loadbal);
551: PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: reuse nonzero structure %sn",jac->ruse ? "true" : "false");
552: PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: print info to screen %sn",jac->logging ? "true" : "false");
553: if (jac->symt == 0) {
554: symt = "nonsymmetric matrix and preconditioner";
555: } else if (jac->symt == 1) {
556: symt = "SPD matrix and preconditioner";
557: } else if (jac->symt == 2) {
558: symt = "nonsymmetric matrix but SPD preconditioner";
559: } else {
560: SETERRQ1(1,"Unknown HYPRE ParaSails sym option %d",jac->symt);
561: }
562: PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: %sn",symt);
563: }
564: return(0);
565: }
566: /* ---------------------------------------------------------------------------------*/
568: static int PCHYPRESetType_HYPRE(PC pc,char *name)
569: {
570: PC_HYPRE *jac = (PC_HYPRE*)pc->data;
571: int ierr;
572: PetscTruth flag;
575: if (pc->ops->setup) {
576: SETERRQ(1,"Cannot set the HYPRE preconditioner type once it has been set");
577: }
579: pc->ops->setup = PCSetUp_HYPRE;
580: pc->ops->apply = PCApply_HYPRE;
581: pc->ops->destroy = PCDestroy_HYPRE;
583: jac->maxiter = PETSC_DEFAULT;
584: jac->tol = PETSC_DEFAULT;
585: jac->printstatistics = PetscLogPrintInfo;
587: PetscStrcmp("pilut",name,&flag);
588: if (flag) {
589: ierr = HYPRE_ParCSRPilutCreate(pc->comm,&jac->hsolver);
590: pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
591: pc->ops->view = PCView_HYPRE_Pilut;
592: jac->destroy = HYPRE_ParCSRPilutDestroy;
593: jac->setup = HYPRE_ParCSRPilutSetup;
594: jac->solve = HYPRE_ParCSRPilutSolve;
595: jac->factorrowsize = PETSC_DEFAULT;
596: return(0);
597: }
598: PetscStrcmp("parasails",name,&flag);
599: if (flag) {
600: ierr = HYPRE_ParaSailsCreate(pc->comm,&jac->hsolver);
601: pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
602: pc->ops->view = PCView_HYPRE_ParaSails;
603: jac->destroy = HYPRE_ParaSailsDestroy;
604: jac->setup = HYPRE_ParaSailsSetup;
605: jac->solve = HYPRE_ParaSailsSolve;
606: return(0);
607: }
608: PetscStrcmp("euclid",name,&flag);
609: if (flag) {
610: ierr = HYPRE_EuclidCreate(pc->comm,&jac->hsolver);
611: pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
612: pc->ops->view = PCView_HYPRE_Euclid;
613: jac->destroy = HYPRE_EuclidDestroy;
614: jac->setup = HYPRE_EuclidSetup;
615: jac->solve = HYPRE_EuclidSolve;
616: return(0);
617: }
618: PetscStrcmp("boomeramg",name,&flag);
619: if (flag) {
620: ierr = HYPRE_BoomerAMGCreate(&jac->hsolver);
621: pc->ops->setfromoptions = PCSetFromOptions_HYPRE_BoomerAMG;
622: pc->ops->view = PCView_HYPRE_BoomerAMG;
623: jac->destroy = HYPRE_BoomerAMGDestroy;
624: jac->setup = HYPRE_BoomerAMGSetup;
625: jac->solve = HYPRE_BoomerAMGSolve;
626: return(0);
627: }
628: SETERRQ1(1,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
629: return(0);
630: }
632: /*
633: It only gets here if the HYPRE type has not been set before the call to
634: ...SetFromOptions() which actually is most of the time
635: */
636: static int PCSetFromOptions_HYPRE(PC pc)
637: {
638: int ierr;
639: char buff[32],*type[] = {"pilut","parasails","boomerAMG","euclid"};
640: PetscTruth flg;
643: PetscOptionsHead("HYPRE preconditioner options");
644: PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"pilut",buff,32,&flg);
646:
647: if (PetscOptionsPublishCount) { /* force the default if it was not yet set and user did not set with option */
648: if (!flg && !pc->ops->apply) {
649: flg = PETSC_TRUE;
650: PetscStrcpy(buff,"pilut");
651: }
652: }
654: if (flg) {
655: PCHYPRESetType_HYPRE(pc,buff);
656: }
657: if (pc->ops->setfromoptions) {
658: pc->ops->setfromoptions(pc);
659: }
660: PetscOptionsTail();
661: return(0);
662: }
664: EXTERN_C_BEGIN
665: int PCCreate_HYPRE(PC pc)
666: {
667: PC_HYPRE *jac;
668: int ierr;
671: ierr = PetscNew(PC_HYPRE,&jac);
672: ierr = PetscMemzero(jac,sizeof(PC_HYPRE));
673: pc->data = jac;
675: pc->ops->setfromoptions = PCSetFromOptions_HYPRE;
676: return(0);
677: }
678: EXTERN_C_END