Actual source code: hypre.c

  1: #define PETSCKSP_DLL

  3: /* 
  4:    Provides an interface to the LLNL package hypre
  5: */

  7: /* Must use hypre 2.0.0 or more recent. */

 9:  #include private/pcimpl.h
 11: #include "HYPRE.h"
 12: #include "HYPRE_parcsr_ls.h"
 13: #include "_hypre_parcsr_mv.h"
 14: #include "_hypre_IJ_mv.h"

 17: EXTERN PetscErrorCode MatHYPRE_IJMatrixCreate(Mat,HYPRE_IJMatrix*);
 18: EXTERN PetscErrorCode MatHYPRE_IJMatrixCopy(Mat,HYPRE_IJMatrix);
 19: EXTERN PetscErrorCode VecHYPRE_IJVectorCreate(Vec,HYPRE_IJVector*);

 21: /* 
 22:    Private context (data structure) for the  preconditioner.  
 23: */
 24: typedef struct {
 25:   HYPRE_Solver       hsolver;
 26:   HYPRE_IJMatrix     ij;
 27:   HYPRE_IJVector     b,x;

 29:   PetscErrorCode     (*destroy)(HYPRE_Solver);
 30:   PetscErrorCode     (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
 31:   PetscErrorCode     (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
 32: 
 33:   MPI_Comm           comm_hypre;
 34:   char              *hypre_type;

 36:   /* options for Pilut and BoomerAMG*/
 37:   int                maxiter;
 38:   double             tol;

 40:   /* options for Pilut */
 41:   int                factorrowsize;

 43:   /* options for ParaSails */
 44:   int                nlevels;
 45:   double             threshhold;
 46:   double             filter;
 47:   int                sym;
 48:   double             loadbal;
 49:   int                logging;
 50:   int                ruse;
 51:   int                symt;

 53:   /* options for Euclid */
 54:   PetscTruth         bjilu;
 55:   int                levels;

 57:   /* options for Euclid and BoomerAMG */
 58:   PetscTruth         printstatistics;

 60:   /* options for BoomerAMG */
 61:   int                cycletype;
 62:   int                maxlevels;
 63:   double             strongthreshold;
 64:   double             maxrowsum;
 65:   int                gridsweeps[3];
 66:   int                coarsentype;
 67:   int                measuretype;
 68:   int                relaxtype[3];
 69:   double             relaxweight;
 70:   double             outerrelaxweight;
 71:   int                relaxorder;
 72:   double             truncfactor;
 73:   PetscTruth         applyrichardson;
 74:   int                pmax;
 75:   int                interptype;
 76:   int                agg_nl;
 77:   int                agg_num_paths;
 78:   int                nodal_coarsen;
 79:   PetscTruth         nodal_relax;
 80:   int                nodal_relax_levels;
 81: } PC_HYPRE;


 86: static PetscErrorCode PCSetUp_HYPRE(PC pc)
 87: {
 88:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
 89:   PetscErrorCode     ierr;
 90:   HYPRE_ParCSRMatrix hmat;
 91:   HYPRE_ParVector    bv,xv;
 92:   PetscInt           bs;
 93:   int                hierr;

 96:   if (!jac->hypre_type) {
 97:     PCHYPRESetType(pc,"pilut");
 98:   }
 99: #if 1
100:   if (!pc->setupcalled) {
101:     /* create the matrix and vectors the first time through */
102:     Vec x,b;
103:     MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);
104:     MatGetVecs(pc->pmat,&x,&b);
105:     VecHYPRE_IJVectorCreate(x,&jac->x);
106:     VecHYPRE_IJVectorCreate(b,&jac->b);
107:     VecDestroy(x);
108:     VecDestroy(b);
109:   } else if (pc->flag != SAME_NONZERO_PATTERN) {
110:     /* rebuild the matrix from scratch */
111:     HYPRE_IJMatrixDestroy(jac->ij);
112:     MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);
113:   }
114: #else  
115:   if (!jac->ij) { /* create the matrix the first time through */
116:     MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);
117:   }
118:   if (!jac->b) { /* create the vectors the first time through */
119:     Vec x,b;
120:     MatGetVecs(pc->pmat,&x,&b);
121:     VecHYPRE_IJVectorCreate(x,&jac->x);
122:     VecHYPRE_IJVectorCreate(b,&jac->b);
123:     VecDestroy(x);
124:     VecDestroy(b);
125:   }
126: #endif
127:   /* special case for BoomerAMG */
128:   if (jac->setup == HYPRE_BoomerAMGSetup) {
129:     MatGetBlockSize(pc->pmat,&bs);
130:     if (bs > 1) {
131:       HYPRE_BoomerAMGSetNumFunctions(jac->hsolver,bs);
132:     }
133:   };
134:   MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);
135:   HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);
136:   HYPRE_IJVectorGetObject(jac->b,(void**)&bv);
137:   HYPRE_IJVectorGetObject(jac->x,(void**)&xv);
138:   h(*jac->setup)(jac->hsolver,hmat,bv,xv);
139:   if (hierr) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE setup, error code %d",hierr);
140:   return(0);
141: }

143: /*
144:     Replaces the address where the HYPRE vector points to its data with the address of
145:   PETSc's data. Saves the old address so it can be reset when we are finished with it.
146:   Allows use to get the data into a HYPRE vector without the cost of memcopies 
147: */
148: #define HYPREReplacePointer(b,newvalue,savedvalue) {\
149:    hypre_ParVector *par_vector   = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));\
150:    hypre_Vector    *local_vector = hypre_ParVectorLocalVector(par_vector);\
151:    savedvalue         = local_vector->data;\
152:    local_vector->data = newvalue;}

156: static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
157: {
158:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
159:   PetscErrorCode     ierr;
160:   HYPRE_ParCSRMatrix hmat;
161:   PetscScalar        *bv,*xv;
162:   HYPRE_ParVector    jbv,jxv;
163:   PetscScalar        *sbv,*sxv;
164:   int                hierr;

167:   if (!jac->applyrichardson) {VecSet(x,0.0);}
168:   VecGetArray(b,&bv);
169:   VecGetArray(x,&xv);
170:   HYPREReplacePointer(jac->b,bv,sbv);
171:   HYPREReplacePointer(jac->x,xv,sxv);

173:   HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);
174:   HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);
175:   HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);
176:   h(*jac->solve)(jac->hsolver,hmat,jbv,jxv);

178:   /*if (hierr && (hierr != HYPRE_ERROR_CONV || jac->solve != HYPRE_BoomerAMGSolve))SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
179:    */
180:  /* error code of HYPRE_ERROR_CONV means convergence not achieved - if
181:     the tolerance is set to 0.0 (the default), a convergence error will
182:     not occur (so we may not want to overide the conv. error here?*/
183:  if (hierr && hierr != HYPRE_ERROR_CONV)
184:   {
185:      SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
186:   }
187:  if (hierr) hypre__global_error = 0;
188: 

190:   HYPREReplacePointer(jac->b,sbv,bv);
191:   HYPREReplacePointer(jac->x,sxv,xv);
192:   VecRestoreArray(x,&xv);
193:   VecRestoreArray(b,&bv);
194:   return(0);
195: }

199: static PetscErrorCode PCDestroy_HYPRE(PC pc)
200: {
201:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;

205:   if (jac->ij) { HYPRE_IJMatrixDestroy(jac->ij); }
206:   if (jac->b)  { HYPRE_IJVectorDestroy(jac->b);  }
207:   if (jac->x)  { HYPRE_IJVectorDestroy(jac->x);  }
208:   if (jac->destroy) { (*jac->destroy)(jac->hsolver); }
209:   PetscStrfree(jac->hypre_type);
210:   if (jac->comm_hypre != MPI_COMM_NULL) { MPI_Comm_free(&(jac->comm_hypre));}
211:   PetscFree(jac);

213:   PetscObjectChangeTypeName((PetscObject)pc,0);
214:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);
215:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);
216:   return(0);
217: }

219: /* --------------------------------------------------------------------------------------------*/
222: static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc)
223: {
224:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
226:   PetscTruth     flag;

229:   PetscOptionsHead("HYPRE Pilut Options");
230:   PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);
231:   if (flag) {
232:     HYPRE_ParCSRPilutSetMaxIter(jac->hsolver,jac->maxiter);
233:   }
234:   PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);
235:   if (flag) {
236:     HYPRE_ParCSRPilutSetDropTolerance(jac->hsolver,jac->tol);
237:   }
238:   PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);
239:   if (flag) {
240:     HYPRE_ParCSRPilutSetFactorRowSize(jac->hsolver,jac->factorrowsize);
241:   }
242:   PetscOptionsTail();
243:   return(0);
244: }

248: static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
249: {
250:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
252:   PetscTruth     iascii;

255:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
256:   if (iascii) {
257:     PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");
258:     if (jac->maxiter != PETSC_DEFAULT) {
259:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);
260:     } else {
261:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");
262:     }
263:     if (jac->tol != PETSC_DEFAULT) {
264:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %G\n",jac->tol);
265:     } else {
266:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");
267:     }
268:     if (jac->factorrowsize != PETSC_DEFAULT) {
269:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);
270:     } else {
271:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");
272:     }
273:   }
274:   return(0);
275: }

277: /* --------------------------------------------------------------------------------------------*/
280: static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc)
281: {
282:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
284:   PetscTruth     flag;
285:   char           *args[8],levels[16];
286:   PetscInt       cnt = 0;

289:   PetscOptionsHead("HYPRE Euclid Options");
290:   PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);
291:   if (flag) {
292:     if (jac->levels < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels);
293:     sprintf(levels,"%d",jac->levels);
294:     args[cnt++] = (char*)"-level"; args[cnt++] = levels;
295:   }
296:   PetscOptionsTruth("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);
297:   if (jac->bjilu) {
298:     args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1";
299:   }
300: 
301:   PetscOptionsTruth("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);
302:   if (jac->printstatistics) {
303:     args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1";
304:     args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1";
305:   }
306:   PetscOptionsTail();
307:   if (cnt) {
308:     HYPRE_EuclidSetParams(jac->hsolver,cnt,args);
309:   }
310:   return(0);
311: }

315: static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
316: {
317:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
319:   PetscTruth     iascii;

322:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
323:   if (iascii) {
324:     PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n");
325:     PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: number of levels %d\n",jac->levels);
326:     if (jac->bjilu) {
327:       PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");
328:     }
329:   }
330:   return(0);
331: }

333: /* --------------------------------------------------------------------------------------------*/

337: static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
338: {
339:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
340:   PetscErrorCode     ierr;
341:   HYPRE_ParCSRMatrix hmat;
342:   PetscScalar        *bv,*xv;
343:   HYPRE_ParVector    jbv,jxv;
344:   PetscScalar        *sbv,*sxv;
345:   int                hierr;

348:   VecSet(x,0.0);
349:   VecGetArray(b,&bv);
350:   VecGetArray(x,&xv);
351:   HYPREReplacePointer(jac->b,bv,sbv);
352:   HYPREReplacePointer(jac->x,xv,sxv);

354:   HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);
355:   HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);
356:   HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);
357: 
358:   hHYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv);
359:   /* error code of 1 in BoomerAMG merely means convergence not achieved */
360:   if (hierr && (hierr != 1)) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
361:   if (hierr) hypre__global_error = 0;
362: 
363:   HYPREReplacePointer(jac->b,sbv,bv);
364:   HYPREReplacePointer(jac->x,sxv,xv);
365:   VecRestoreArray(x,&xv);
366:   VecRestoreArray(b,&bv);
367:   return(0);
368: }

370: static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
371: static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"};
372: static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
373: static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","","SOR/Jacobi","backward-SOR/Jacobi","","symmetric-SOR/Jacobi",
374:                                                   "","","Gaussian-elimination"};
375: static const char *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i",
376:                                                   "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"};
379: static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
380: {
381:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
383:   int            n,indx;
384:   PetscTruth     flg, tmp_truth;
385:   double         tmpdbl, twodbl[2];

388:   PetscOptionsHead("HYPRE BoomerAMG Options");
389:   PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);
390:   if (flg) {
391:     jac->cycletype = indx;
392:     HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);
393:   }
394:   PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);
395:   if (flg) {
396:     if (jac->maxlevels < 2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
397:     HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);
398:   }
399:   PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);
400:   if (flg) {
401:     if (jac->maxiter < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
402:     HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);
403:   }
404:   PetscOptionsScalar("-pc_hypre_boomeramg_tol","Convergence tolerance PER hypre call (0.0 = use a fixed number of iterations)","None",jac->tol,&jac->tol,&flg);
405:   if (flg) {
406:     if (jac->tol < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be greater than or equal to zero",jac->tol);
407:     HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);
408:   }

410:   PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);
411:   if (flg) {
412:     if (jac->truncfactor < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor);
413:     HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);
414:   }

416:  PetscOptionsInt("-pc_hypre_boomeramg_P_max","Max elements per row for interpolation operator ( 0=unlimited )","None",jac->pmax,&jac->pmax,&flg);
417:   if (flg) {
418:     if (jac->pmax < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"P_max %G must be greater than or equal to zero",jac->pmax);
419:     HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);
420:   }

422:  PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);
423:   if (flg) {
424:      if (jac->agg_nl < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %G must be greater than or equal to zero",jac->agg_nl);

426:      HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl);
427:   }


430:  PetscOptionsInt("-pc_hypre_boomeramg_agg_num_paths","Number of paths for aggressive coarsening","None",jac->agg_num_paths,&jac->agg_num_paths,&flg);
431:   if (flg) {
432:      if (jac->agg_num_paths < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of paths %G must be greater than or equal to 1",jac->agg_num_paths);

434:      HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);
435:   }


438:   PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);
439:   if (flg) {
440:     if (jac->strongthreshold < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold);
441:     HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);
442:   }
443:   PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);
444:   if (flg) {
445:     if (jac->maxrowsum < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum);
446:     if (jac->maxrowsum > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum);
447:     HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);
448:   }

450:   /* Grid sweeps */
451:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_all","Number of sweeps for the up and down grid levels","None",jac->gridsweeps[0],&indx,&flg);
452:   if (flg) {
453:     HYPRE_BoomerAMGSetNumSweeps(jac->hsolver,indx);
454:     /* modify the jac structure so we can view the updated options with PC_View */
455:     jac->gridsweeps[0] = indx;
456:     jac->gridsweeps[1] = indx;
457:     /*defaults coarse to 1 */
458:     jac->gridsweeps[2] = 1;
459:   }

461:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx ,&flg);
462:   if (flg) {
463:     HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 1);
464:     jac->gridsweeps[0] = indx;
465:   }
466:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);
467:   if (flg) {
468:     HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 2);
469:     jac->gridsweeps[1] = indx;
470:   }
471:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);
472:   if (flg) {
473:     HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 3);
474:     jac->gridsweeps[2] = indx;
475:   }

477:   /* Relax type */
478:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_all","Relax type for the up and down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);
479:   if (flg) {
480:     jac->relaxtype[0] = jac->relaxtype[1]  = indx;
481:     HYPRE_BoomerAMGSetRelaxType(jac->hsolver, indx);
482:     /* by default, coarse type set to 9 */
483:     jac->relaxtype[2] = 9;
484: 
485:   }
486:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);
487:   if (flg) {
488:     jac->relaxtype[0] = indx;
489:     HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 1);
490:   }
491:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);
492:   if (flg) {
493:     jac->relaxtype[1] = indx;
494:     HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 2);
495:   }
496:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);
497:   if (flg) {
498:     jac->relaxtype[2] = indx;
499:     HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 3);
500:   }

502:   /* Relaxation Weight */
503:   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);
504:   if (flg) {
505:     HYPRE_BoomerAMGSetRelaxWt(jac->hsolver,tmpdbl);
506:     jac->relaxweight = tmpdbl;
507:   }

509:   n=2;
510:   twodbl[0] = twodbl[1] = 1.0;
511:   PetscOptionsRealArray("-pc_hypre_boomeramg_relax_weight_level","Set the relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);
512:   if (flg) {
513:     if (n == 2) {
514:       indx =  (int)PetscAbsReal(twodbl[1]);
515:       HYPRE_BoomerAMGSetLevelRelaxWt(jac->hsolver,twodbl[0],indx);
516:     } else {
517:       SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Relax weight level: you must provide 2 values separated by a comma (and no space), you provided %d",n);
518:     }
519:   }

521:   /* Outer relaxation Weight */
522:   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);
523:   if (flg) {
524:     HYPRE_BoomerAMGSetOuterWt( jac->hsolver, tmpdbl);
525:     jac->outerrelaxweight = tmpdbl;
526:   }

528:   n=2;
529:   twodbl[0] = twodbl[1] = 1.0;
530:   PetscOptionsRealArray("-pc_hypre_boomeramg_outer_relax_weight_level","Set the outer relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);
531:   if (flg) {
532:     if (n == 2) {
533:       indx =  (int)PetscAbsReal(twodbl[1]);
534:       HYPRE_BoomerAMGSetLevelOuterWt( jac->hsolver, twodbl[0], indx);
535:     } else {
536:       SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Relax weight outer level: You must provide 2 values separated by a comma (and no space), you provided %d",n);
537:     }
538:   }

540:   /* the Relax Order */
541:   /* PetscOptionsName("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", &flg); */
542:   PetscOptionsTruth( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);

544:   if (flg) {
545:     jac->relaxorder = 0;
546:     HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);
547:   }
548:   PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);
549:   if (flg) {
550:     jac->measuretype = indx;
551:     HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);
552:   }
553:   /* update list length 3/07 */
554:   PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,11,HYPREBoomerAMGCoarsenType[6],&indx,&flg);
555:   if (flg) {
556:     jac->coarsentype = indx;
557:     HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);
558:   }
559: 
560:   /* new 3/07 */
561:   PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,14,HYPREBoomerAMGInterpType[0],&indx,&flg);
562:   if (flg) {
563:     jac->interptype = indx;
564:     HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);
565:   }



569:   PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);
570:   if (flg) {
571:     int level=3;
572:     jac->printstatistics = PETSC_TRUE;
573:     PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);
574:     HYPRE_BoomerAMGSetPrintLevel(jac->hsolver,level);
575:     HYPRE_BoomerAMGSetDebugFlag(jac->hsolver,level);
576:   }

578:   PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", PETSC_FALSE, &tmp_truth, &flg);
579:   if (flg && tmp_truth) {
580:     jac->nodal_coarsen = 1;
581:     HYPRE_BoomerAMGSetNodal(jac->hsolver,1);
582:   }

584:   PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);
585:   if (flg && tmp_truth) {
586:     PetscInt tmp_int;
587:     PetscOptionsInt( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);
588:     if (flg) jac->nodal_relax_levels = tmp_int;
589:     HYPRE_BoomerAMGSetSmoothType(jac->hsolver,6);
590:     HYPRE_BoomerAMGSetDomainType(jac->hsolver,1);
591:     HYPRE_BoomerAMGSetOverlap(jac->hsolver,0);
592:     HYPRE_BoomerAMGSetSmoothNumLevels(jac->hsolver,jac->nodal_relax_levels);
593:   }

595:   PetscOptionsTail();
596:   return(0);
597: }

601: static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
602: {
603:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;

607:   HYPRE_BoomerAMGSetMaxIter(jac->hsolver,its*jac->maxiter);
608:   HYPRE_BoomerAMGSetTol(jac->hsolver,PetscMin(rtol,jac->tol));
609:   jac->applyrichardson = PETSC_TRUE;
610:   PCApply_HYPRE(pc,b,y);
611:   jac->applyrichardson = PETSC_FALSE;
612:   HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);
613:   HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);
614:   return(0);
615: }


620: static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
621: {
622:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
624:   PetscTruth     iascii;

627:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
628:   if (iascii) {
629:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");
630:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);
631:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);
632:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);
633:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);
634:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);
635:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation truncation factor %G\n",jac->truncfactor);
636:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);
637:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);
638:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);
639: 
640:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);

642:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps down         %d\n",jac->gridsweeps[0]);
643:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps up           %d\n",jac->gridsweeps[1]);
644:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on coarse    %d\n",jac->gridsweeps[2]);

646:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);
647:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);
648:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);

650:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax weight  (all)      %G\n",jac->relaxweight);
651:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);

653:     if (jac->relaxorder) {
654:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using CF-relaxation\n");
655:     } else {
656:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using CF-relaxation\n");
657:     }
658:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);
659:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);
660:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation type  %s\n",HYPREBoomerAMGInterpType[jac->interptype]);
661:     if (jac->nodal_coarsen) {
662:       PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal())\n");
663:     }
664:     if (jac->nodal_relax) {
665:       PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);
666:     }
667:   }
668:   return(0);
669: }

671: /* --------------------------------------------------------------------------------------------*/
674: static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc)
675: {
676:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
678:   int            indx;
679:   PetscTruth     flag;
680:   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};

683:   PetscOptionsHead("HYPRE ParaSails Options");
684:   PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);
685:   PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);
686:   if (flag) {
687:     HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);
688:   }

690:   PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);
691:   if (flag) {
692:     HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);
693:   }

695:   PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);
696:   if (flag) {
697:     HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);
698:   }

700:   PetscOptionsTruth("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,&flag);
701:   if (flag) {
702:     HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);
703:   }

705:   PetscOptionsTruth("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,&flag);
706:   if (flag) {
707:     HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);
708:   }

710:   PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);
711:   if (flag) {
712:     jac->symt = indx;
713:     HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);
714:   }

716:   PetscOptionsTail();
717:   return(0);
718: }

722: static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
723: {
724:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
726:   PetscTruth     iascii;
727:   const char     *symt = 0;;

730:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
731:   if (iascii) {
732:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");
733:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);
734:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %G\n",jac->threshhold);
735:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %G\n",jac->filter);
736:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %G\n",jac->loadbal);
737:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscTruths[jac->ruse]);
738:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscTruths[jac->logging]);
739:     if (!jac->symt) {
740:       symt = "nonsymmetric matrix and preconditioner";
741:     } else if (jac->symt == 1) {
742:       symt = "SPD matrix and preconditioner";
743:     } else if (jac->symt == 2) {
744:       symt = "nonsymmetric matrix but SPD preconditioner";
745:     } else {
746:       SETERRQ1(PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
747:     }
748:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);
749:   }
750:   return(0);
751: }
752: /* ---------------------------------------------------------------------------------*/

757: PetscErrorCode  PCHYPREGetType_HYPRE(PC pc,const char *name[])
758: {
759:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;

762:   *name = jac->hypre_type;
763:   return(0);
764: }

770: PetscErrorCode  PCHYPRESetType_HYPRE(PC pc,const char name[])
771: {
772:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
774:   PetscTruth     flag;

777:   if (jac->hypre_type) {
778:     PetscStrcmp(jac->hypre_type,name,&flag);
779:     if (!flag) SETERRQ(PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
780:     return(0);
781:   } else {
782:     PetscStrallocpy(name, &jac->hypre_type);
783:   }
784: 
785:   jac->maxiter            = PETSC_DEFAULT;
786:   jac->tol                = PETSC_DEFAULT;
787:   jac->printstatistics    = PetscLogPrintInfo;

789:   PetscStrcmp("pilut",jac->hypre_type,&flag);
790:   if (flag) {
791:     HYPRE_ParCSRPilutCreate(jac->comm_hypre,&jac->hsolver);
792:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
793:     pc->ops->view           = PCView_HYPRE_Pilut;
794:     jac->destroy            = HYPRE_ParCSRPilutDestroy;
795:     jac->setup              = HYPRE_ParCSRPilutSetup;
796:     jac->solve              = HYPRE_ParCSRPilutSolve;
797:     jac->factorrowsize      = PETSC_DEFAULT;
798:     return(0);
799:   }
800:   PetscStrcmp("parasails",jac->hypre_type,&flag);
801:   if (flag) {
802:     HYPRE_ParaSailsCreate(jac->comm_hypre,&jac->hsolver);
803:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
804:     pc->ops->view           = PCView_HYPRE_ParaSails;
805:     jac->destroy            = HYPRE_ParaSailsDestroy;
806:     jac->setup              = HYPRE_ParaSailsSetup;
807:     jac->solve              = HYPRE_ParaSailsSolve;
808:     /* initialize */
809:     jac->nlevels     = 1;
810:     jac->threshhold  = .1;
811:     jac->filter      = .1;
812:     jac->loadbal     = 0;
813:     if (PetscLogPrintInfo) {
814:       jac->logging     = (int) PETSC_TRUE;
815:     } else {
816:       jac->logging     = (int) PETSC_FALSE;
817:     }
818:     jac->ruse = (int) PETSC_FALSE;
819:     jac->symt   = 0;
820:     HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);
821:     HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);
822:     HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);
823:     HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);
824:     HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);
825:     HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);
826:     return(0);
827:   }
828:   PetscStrcmp("euclid",jac->hypre_type,&flag);
829:   if (flag) {
830:     HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver);
831:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
832:     pc->ops->view           = PCView_HYPRE_Euclid;
833:     jac->destroy            = HYPRE_EuclidDestroy;
834:     jac->setup              = HYPRE_EuclidSetup;
835:     jac->solve              = HYPRE_EuclidSolve;
836:     /* initialization */
837:     jac->bjilu              = PETSC_FALSE;
838:     jac->levels             = 1;
839:     return(0);
840:   }
841:   PetscStrcmp("boomeramg",jac->hypre_type,&flag);
842:   if (flag) {
843:     HYPRE_BoomerAMGCreate(&jac->hsolver);
844:     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
845:     pc->ops->view            = PCView_HYPRE_BoomerAMG;
846:     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
847:     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
848:     jac->destroy             = HYPRE_BoomerAMGDestroy;
849:     jac->setup               = HYPRE_BoomerAMGSetup;
850:     jac->solve               = HYPRE_BoomerAMGSolve;
851:     jac->applyrichardson     = PETSC_FALSE;
852:     /* these defaults match the hypre defaults */
853:     jac->cycletype        = 1;
854:     jac->maxlevels        = 25;
855:     jac->maxiter          = 1;
856:     jac->tol              = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */
857:     jac->truncfactor      = 0.0;
858:     jac->strongthreshold  = .25;
859:     jac->maxrowsum        = .9;
860:     jac->coarsentype      = 6;
861:     jac->measuretype      = 0;
862:     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = 1;
863:     jac->relaxtype[0]     = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */
864:     jac->relaxtype[2]     = 9; /*G.E. */
865:     jac->relaxweight      = 1.0;
866:     jac->outerrelaxweight = 1.0;
867:     jac->relaxorder       = 1;
868:     jac->interptype       = 0;
869:     jac->agg_nl           = 0;
870:     jac->pmax             = 0;
871:     jac->truncfactor      = 0.0;
872:     jac->agg_num_paths    = 1;

874:     jac->nodal_coarsen    = 0;
875:     jac->nodal_relax      = PETSC_FALSE;
876:     jac->nodal_relax_levels = 1;
877:     HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);
878:     HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);
879:     HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);
880:     HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);
881:     HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);
882:     HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);
883:     HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);
884:     HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);
885:     HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);
886:     HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);
887:     HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);
888:     HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl);
889:     HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);
890:     HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);
891:     HYPRE_BoomerAMGSetRelaxType(jac->hsolver, jac->relaxtype[0]);  /*defaults coarse to 9*/
892:     HYPRE_BoomerAMGSetNumSweeps(jac->hsolver, jac->gridsweeps[0]); /*defaults coarse to 1 */
893:     return(0);
894:   }
895:   PetscStrfree(jac->hypre_type);
896:   jac->hypre_type = PETSC_NULL;
897:   SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
898:   return(0);
899: }

902: /*
903:     It only gets here if the HYPRE type has not been set before the call to 
904:    ...SetFromOptions() which actually is most of the time
905: */
908: static PetscErrorCode PCSetFromOptions_HYPRE(PC pc)
909: {
910:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
912:   int            indx;
913:   const char     *type[] = {"pilut","parasails","boomeramg","euclid"};
914:   PetscTruth     flg;

917:   PetscOptionsHead("HYPRE preconditioner options");
918:   PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"pilut",&indx,&flg);
919:   if (PetscOptionsPublishCount) {   /* force the default if it was not yet set and user did not set with option */
920:     if (!flg && !jac->hypre_type) {
921:       flg   = PETSC_TRUE;
922:       indx = 0;
923:     }
924:   }
925:   if (flg) {
926:     PCHYPRESetType_HYPRE(pc,type[indx]);
927:   }
928:   if (pc->ops->setfromoptions) {
929:     pc->ops->setfromoptions(pc);
930:   }
931:   PetscOptionsTail();
932:   return(0);
933: }

937: /*@C
938:      PCHYPRESetType - Sets which hypre preconditioner you wish to use

940:    Input Parameters:
941: +     pc - the preconditioner context
942: -     name - either  pilut, parasails, boomeramg, euclid

944:    Options Database Keys:
945:    -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
946:  
947:    Level: intermediate

949: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
950:            PCHYPRE

952: @*/
953: PetscErrorCode  PCHYPRESetType(PC pc,const char name[])
954: {
955:   PetscErrorCode ierr,(*f)(PC,const char[]);

960:   PetscObjectQueryFunction((PetscObject)pc,"PCHYPRESetType_C",(void (**)(void))&f);
961:   if (f) {
962:     (*f)(pc,name);
963:   }
964:   return(0);
965: }

969: /*@C
970:      PCHYPREGetType - Gets which hypre preconditioner you are using

972:    Input Parameter:
973: .     pc - the preconditioner context

975:    Output Parameter:
976: .     name - either  pilut, parasails, boomeramg, euclid

978:    Level: intermediate

980: .seealso:  PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC,
981:            PCHYPRE

983: @*/
984: PetscErrorCode  PCHYPREGetType(PC pc,const char *name[])
985: {
986:   PetscErrorCode ierr,(*f)(PC,const char*[]);

991:   PetscObjectQueryFunction((PetscObject)pc,"PCHYPREGetType_C",(void (**)(void))&f);
992:   if (f) {
993:     (*f)(pc,name);
994:   }
995:   return(0);
996: }

998: /*MC
999:      PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre

1001:    Options Database Keys:
1002: +   -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
1003: -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
1004:           preconditioner
1005:  
1006:    Level: intermediate

1008:    Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()),
1009:           the many hypre options can ONLY be set via the options database (e.g. the command line
1010:           or with PetscOptionsSetValue(), there are no functions to set them)

1012:           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations
1013:           (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if 
1014:           -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner 
1015:           (-pc_hypre_boomeramg_rtol should be set to 0.0 - the default - to strictly use a fixed number of 
1016:           iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations 
1017:           and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10 
1018:           then AT MOST twenty V-cycles of boomeramg will be called.

1020:            Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation 
1021:            (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry.  
1022:            Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi.
1023:           If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly
1024:           and use -ksp_max_it to control the number of V-cycles.
1025:           (see the PETSc FAQ.html at the PETSc website under the Documentation tab).

1027:           2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option
1028:           -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L.

1030: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1031:            PCHYPRESetType()

1033: M*/

1038: PetscErrorCode  PCCreate_HYPRE(PC pc)
1039: {
1040:   PC_HYPRE       *jac;

1044:   PetscNew(PC_HYPRE,&jac);
1045:   PetscLogObjectMemory(pc,sizeof(PC_HYPRE));
1046:   pc->data                 = jac;
1047:   pc->ops->destroy         = PCDestroy_HYPRE;
1048:   pc->ops->setfromoptions  = PCSetFromOptions_HYPRE;
1049:   pc->ops->setup           = PCSetUp_HYPRE;
1050:   pc->ops->apply           = PCApply_HYPRE;
1051:   jac->comm_hypre          = MPI_COMM_NULL;
1052:   jac->hypre_type          = PETSC_NULL;
1053:   /* duplicate communicator for hypre */
1054:   MPI_Comm_dup(pc->comm,&(jac->comm_hypre));
1055:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);
1056:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);
1057:   return(0);
1058: }