Actual source code: hypre.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  2: /*
  3:    Provides an interface to the LLNL package hypre
  4: */

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

  8: #include <petsc/private/pcimpl.h>          /*I "petscpc.h" I*/
  9: #include <../src/dm/impls/da/hypre/mhyp.h>
 10: #include <_hypre_parcsr_ls.h>

 12: static PetscBool cite = PETSC_FALSE;
 13: static const char hypreCitation[] = "@manual{hypre-web-page,\n  title  = {{\\sl hypre}: High Performance Preconditioners},\n  organization = {Lawrence Livermore National Laboratory},\n  note  = {\\url{http://www.llnl.gov/CASC/hypre/}}\n}\n";

 15: /*
 16:    Private context (data structure) for the  preconditioner.
 17: */
 18: typedef struct {
 19:   HYPRE_Solver   hsolver;
 20:   HYPRE_IJMatrix ij;
 21:   HYPRE_IJVector b,x;

 23:   HYPRE_Int (*destroy)(HYPRE_Solver);
 24:   HYPRE_Int (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
 25:   HYPRE_Int (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
 26:   HYPRE_Int (*setdgrad)(HYPRE_Solver,HYPRE_ParCSRMatrix);
 27:   HYPRE_Int (*setdcurl)(HYPRE_Solver,HYPRE_ParCSRMatrix);
 28:   HYPRE_Int (*setcoord)(HYPRE_Solver,HYPRE_ParVector,HYPRE_ParVector,HYPRE_ParVector);
 29:   HYPRE_Int (*setdim)(HYPRE_Solver,HYPRE_Int);

 31:   MPI_Comm comm_hypre;
 32:   char     *hypre_type;

 34:   /* options for Pilut and BoomerAMG*/
 35:   PetscInt maxiter;
 36:   double   tol;

 38:   /* options for Pilut */
 39:   PetscInt factorrowsize;

 41:   /* options for ParaSails */
 42:   PetscInt nlevels;
 43:   double   threshhold;
 44:   double   filter;
 45:   PetscInt sym;
 46:   double   loadbal;
 47:   PetscInt logging;
 48:   PetscInt ruse;
 49:   PetscInt symt;

 51:   /* options for BoomerAMG */
 52:   PetscBool printstatistics;

 54:   /* options for BoomerAMG */
 55:   PetscInt  cycletype;
 56:   PetscInt  maxlevels;
 57:   double    strongthreshold;
 58:   double    maxrowsum;
 59:   PetscInt  gridsweeps[3];
 60:   PetscInt  coarsentype;
 61:   PetscInt  measuretype;
 62:   PetscInt  relaxtype[3];
 63:   double    relaxweight;
 64:   double    outerrelaxweight;
 65:   PetscInt  relaxorder;
 66:   double    truncfactor;
 67:   PetscBool applyrichardson;
 68:   PetscInt  pmax;
 69:   PetscInt  interptype;
 70:   PetscInt  agg_nl;
 71:   PetscInt  agg_num_paths;
 72:   PetscInt  nodal_coarsen;
 73:   PetscBool nodal_relax;
 74:   PetscInt  nodal_relax_levels;

 76:   /* options for AS (Auxiliary Space preconditioners) */
 77:   PetscInt  as_print;
 78:   PetscInt  as_max_iter;
 79:   PetscReal as_tol;
 80:   PetscInt  as_relax_type;
 81:   PetscInt  as_relax_times;
 82:   PetscReal as_relax_weight;
 83:   PetscReal as_omega;
 84:   PetscInt  as_amg_alpha_opts[5]; /* AMG coarsen type, agg_levels, relax_type, interp_type, Pmax for vector Poisson (AMS) or Curl problem (ADS) */
 85:   PetscReal as_amg_alpha_theta;   /* AMG strength for vector Poisson (AMS) or Curl problem (ADS) */
 86:   PetscInt  as_amg_beta_opts[5];  /* AMG coarsen type, agg_levels, relax_type, interp_type, Pmax for scalar Poisson (AMS) or vector Poisson (ADS) */
 87:   PetscReal as_amg_beta_theta;    /* AMG strength for scalar Poisson (AMS) or vector Poisson (ADS)  */
 88:   PetscInt  ams_cycle_type;
 89:   PetscInt  ads_cycle_type;

 91:   /* additional data */
 92:   HYPRE_IJVector coords[3];
 93:   HYPRE_IJVector constants[3];
 94:   HYPRE_IJMatrix G;
 95:   HYPRE_IJMatrix C;
 96:   HYPRE_IJMatrix alpha_Poisson;
 97:   HYPRE_IJMatrix beta_Poisson;
 98:   PetscBool      ams_beta_is_zero;
 99: } PC_HYPRE;

103: PetscErrorCode PCHYPREGetSolver(PC pc,HYPRE_Solver *hsolver)
104: {
105:   PC_HYPRE *jac = (PC_HYPRE*)pc->data;

108:   *hsolver = jac->hsolver;
109:   return(0);
110: }

114: static PetscErrorCode PCSetUp_HYPRE(PC pc)
115: {
116:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
117:   PetscErrorCode     ierr;
118:   HYPRE_ParCSRMatrix hmat;
119:   HYPRE_ParVector    bv,xv;
120:   PetscInt           bs;

123:   if (!jac->hypre_type) {
124:     PCHYPRESetType(pc,"boomeramg");
125:   }

127:   if (pc->setupcalled) {
128:     /* always destroy the old matrix and create a new memory;
129:        hope this does not churn the memory too much. The problem
130:        is I do not know if it is possible to put the matrix back to
131:        its initial state so that we can directly copy the values
132:        the second time through. */
133:     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->ij));
134:     jac->ij = 0;
135:   }

137:   if (!jac->ij) { /* create the matrix the first time through */
138:     MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);
139:   }
140:   if (!jac->b) { /* create the vectors the first time through */
141:     Vec x,b;
142:     MatCreateVecs(pc->pmat,&x,&b);
143:     VecHYPRE_IJVectorCreate(x,&jac->x);
144:     VecHYPRE_IJVectorCreate(b,&jac->b);
145:     VecDestroy(&x);
146:     VecDestroy(&b);
147:   }

149:   /* special case for BoomerAMG */
150:   if (jac->setup == HYPRE_BoomerAMGSetup) {
151:     MatGetBlockSize(pc->pmat,&bs);
152:     if (bs > 1) PetscStackCallStandard(HYPRE_BoomerAMGSetNumFunctions,(jac->hsolver,bs));
153:   }

155:   /* special case for AMS */
156:   if (jac->setup == HYPRE_AMSSetup) {
157:     if (!jac->coords[0] && !jac->constants[0]) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE AMS preconditioner needs either coordinate vectors via PCSetCoordinates() or edge constant vectors via PCHYPRESetEdgeConstantVectors()");
158:     if (!jac->G) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE AMS preconditioner needs discrete gradient operator via PCHYPRESetDiscreteGradient");
159:   }
160:   /* special case for ADS */
161:   if (jac->setup == HYPRE_ADSSetup) {
162:     if (!jac->coords[0]) {
163:       SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE ADS preconditioner needs coordinate vectors via PCSetCoordinates()");
164:     } else if (!jac->coords[1] || !jac->coords[2]) {
165:       SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE ADS preconditioner has been designed for three dimensional problems! For two dimensional problems, use HYPRE AMS instead");
166:     }
167:     if (!jac->G) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE ADS preconditioner needs discrete gradient operator via PCHYPRESetDiscreteGradient");
168:     if (!jac->C) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE ADS preconditioner needs discrete curl operator via PCHYPRESetDiscreteGradient");
169:   }
170:   MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);
171:   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
172:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&bv));
173:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&xv));
174:   PetscStackCall("HYPRE_SetupXXX",(*jac->setup)(jac->hsolver,hmat,bv,xv););
175:   return(0);
176: }

178: /*
179:     Replaces the address where the HYPRE vector points to its data with the address of
180:   PETSc's data. Saves the old address so it can be reset when we are finished with it.
181:   Allows use to get the data into a HYPRE vector without the cost of memcopies
182: */
183: #define HYPREReplacePointer(b,newvalue,savedvalue) { \
184:     hypre_ParVector *par_vector   = (hypre_ParVector*)hypre_IJVectorObject(((hypre_IJVector*)b)); \
185:     hypre_Vector    *local_vector = hypre_ParVectorLocalVector(par_vector); \
186:     savedvalue         = local_vector->data; \
187:     local_vector->data = newvalue;          \
188: }

192: static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
193: {
194:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
195:   PetscErrorCode     ierr;
196:   HYPRE_ParCSRMatrix hmat;
197:   PetscScalar        *xv;
198:   const PetscScalar  *bv,*sbv;
199:   HYPRE_ParVector    jbv,jxv;
200:   PetscScalar        *sxv;
201:   PetscInt           hierr;

204:   PetscCitationsRegister(hypreCitation,&cite);
205:   if (!jac->applyrichardson) {VecSet(x,0.0);}
206:   VecGetArrayRead(b,&bv);
207:   VecGetArray(x,&xv);
208:   HYPREReplacePointer(jac->b,(PetscScalar*)bv,sbv);
209:   HYPREReplacePointer(jac->x,xv,sxv);

211:   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
212:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&jbv));
213:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&jxv));
214:   PetscStackCall("Hypre solve",h(*jac->solve)(jac->hsolver,hmat,jbv,jxv);
215:                                if (hierr && hierr != HYPRE_ERROR_CONV) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
216:                                if (hierr) hypre__global_error = 0;);

218:   HYPREReplacePointer(jac->b,(PetscScalar*)sbv,bv);
219:   HYPREReplacePointer(jac->x,sxv,xv);
220:   VecRestoreArray(x,&xv);
221:   VecRestoreArrayRead(b,&bv);
222:   return(0);
223: }

227: static PetscErrorCode PCDestroy_HYPRE(PC pc)
228: {
229:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;

233:   if (jac->ij) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->ij));
234:   if (jac->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->b));
235:   if (jac->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->x));
236:   if (jac->coords[0]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[0]));
237:   if (jac->coords[1]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[1]));
238:   if (jac->coords[2]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[2]));
239:   if (jac->constants[0]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[0]));
240:   if (jac->constants[1]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[1]));
241:   if (jac->constants[2]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[2]));
242:   if (jac->G) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->G));
243:   if (jac->C) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->C));
244:   if (jac->alpha_Poisson) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->alpha_Poisson));
245:   if (jac->beta_Poisson) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->beta_Poisson));
246:   if (jac->destroy) PetscStackCall("HYPRE_DestroyXXX",(*jac->destroy)(jac->hsolver););
247:   PetscFree(jac->hypre_type);
248:   if (jac->comm_hypre != MPI_COMM_NULL) { MPI_Comm_free(&(jac->comm_hypre));}
249:   PetscFree(pc->data);

251:   PetscObjectChangeTypeName((PetscObject)pc,0);
252:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetType_C",NULL);
253:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPREGetType_C",NULL);
254:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetCoordinates_C",NULL);
255:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteGradient_C",NULL);
256:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteCurl_C",NULL);
257:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetConstantEdgeVectors_C",NULL);
258:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetAlphaPoissonMatrix_C",NULL);
259:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetBetaPoissonMatrix_C",NULL);
260:   return(0);
261: }

263: /* --------------------------------------------------------------------------------------------*/
266: static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PetscOptions *PetscOptionsObject,PC pc)
267: {
268:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
270:   PetscBool      flag;

273:   PetscOptionsHead(PetscOptionsObject,"HYPRE Pilut Options");
274:   PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);
275:   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetMaxIter,(jac->hsolver,jac->maxiter));
276:   PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);
277:   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetDropTolerance,(jac->hsolver,jac->tol));
278:   PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);
279:   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetFactorRowSize,(jac->hsolver,jac->factorrowsize));
280:   PetscOptionsTail();
281:   return(0);
282: }

286: static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
287: {
288:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
290:   PetscBool      iascii;

293:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
294:   if (iascii) {
295:     PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");
296:     if (jac->maxiter != PETSC_DEFAULT) {
297:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);
298:     } else {
299:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");
300:     }
301:     if (jac->tol != PETSC_DEFAULT) {
302:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %g\n",(double)jac->tol);
303:     } else {
304:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");
305:     }
306:     if (jac->factorrowsize != PETSC_DEFAULT) {
307:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);
308:     } else {
309:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");
310:     }
311:   }
312:   return(0);
313: }

315: /* --------------------------------------------------------------------------------------------*/

319: static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
320: {
321:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
322:   PetscErrorCode     ierr;
323:   HYPRE_ParCSRMatrix hmat;
324:   PetscScalar        *xv;
325:   const PetscScalar  *bv;
326:   HYPRE_ParVector    jbv,jxv;
327:   PetscScalar        *sbv,*sxv;
328:   PetscInt           hierr;

331:   PetscCitationsRegister(hypreCitation,&cite);
332:   VecSet(x,0.0);
333:   VecGetArrayRead(b,&bv);
334:   VecGetArray(x,&xv);
335:   HYPREReplacePointer(jac->b,(PetscScalar*)bv,sbv);
336:   HYPREReplacePointer(jac->x,xv,sxv);

338:   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
339:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&jbv));
340:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&jxv));

342:   hHYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv);
343:   /* error code of 1 in BoomerAMG merely means convergence not achieved */
344:   if (hierr && (hierr != 1)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
345:   if (hierr) hypre__global_error = 0;

347:   HYPREReplacePointer(jac->b,sbv,bv);
348:   HYPREReplacePointer(jac->x,sxv,xv);
349:   VecRestoreArray(x,&xv);
350:   VecRestoreArrayRead(b,&bv);
351:   return(0);
352: }

354: /* static array length */
355: #define ALEN(a) (sizeof(a)/sizeof((a)[0]))

357: static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
358: static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"};
359: static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
360: /* The following corresponds to HYPRE_BoomerAMGSetRelaxType which has many missing numbers in the enum */
361: static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","seqboundary-Gauss-Seidel","SOR/Jacobi","backward-SOR/Jacobi",
362:                                                   "" /* [5] hybrid chaotic Gauss-Seidel (works only with OpenMP) */,"symmetric-SOR/Jacobi",
363:                                                   "" /* 7 */,"l1scaled-SOR/Jacobi","Gaussian-elimination",
364:                                                   "" /* 10 */, "" /* 11 */, "" /* 12 */, "" /* 13 */, "" /* 14 */,
365:                                                   "CG" /* non-stationary */,"Chebyshev","FCF-Jacobi","l1scaled-Jacobi"};
366: static const char *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i",
367:                                                   "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"};
370: static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PetscOptions *PetscOptionsObject,PC pc)
371: {
372:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
374:   PetscInt       n,indx,level;
375:   PetscBool      flg, tmp_truth;
376:   double         tmpdbl, twodbl[2];

379:   PetscOptionsHead(PetscOptionsObject,"HYPRE BoomerAMG Options");
380:   PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType+1,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);
381:   if (flg) {
382:     jac->cycletype = indx+1;
383:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
384:   }
385:   PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);
386:   if (flg) {
387:     if (jac->maxlevels < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
388:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
389:   }
390:   PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);
391:   if (flg) {
392:     if (jac->maxiter < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
393:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
394:   }
395:   PetscOptionsScalar("-pc_hypre_boomeramg_tol","Convergence tolerance PER hypre call (0.0 = use a fixed number of iterations)","None",jac->tol,&jac->tol,&flg);
396:   if (flg) {
397:     if (jac->tol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %g must be greater than or equal to zero",(double)jac->tol);
398:     PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
399:   }

401:   PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);
402:   if (flg) {
403:     if (jac->truncfactor < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %g must be great than or equal zero",(double)jac->truncfactor);
404:     PetscStackCallStandard(HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
405:   }

407:   PetscOptionsInt("-pc_hypre_boomeramg_P_max","Max elements per row for interpolation operator (0=unlimited)","None",jac->pmax,&jac->pmax,&flg);
408:   if (flg) {
409:     if (jac->pmax < 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"P_max %g must be greater than or equal to zero",(double)jac->pmax);
410:     PetscStackCallStandard(HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
411:   }

413:   PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);
414:   if (flg) {
415:     if (jac->agg_nl < 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %g must be greater than or equal to zero",(double)jac->agg_nl);

417:     PetscStackCallStandard(HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
418:   }


421:   PetscOptionsInt("-pc_hypre_boomeramg_agg_num_paths","Number of paths for aggressive coarsening","None",jac->agg_num_paths,&jac->agg_num_paths,&flg);
422:   if (flg) {
423:     if (jac->agg_num_paths < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of paths %g must be greater than or equal to 1",(double)jac->agg_num_paths);

425:     PetscStackCallStandard(HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
426:   }


429:   PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);
430:   if (flg) {
431:     if (jac->strongthreshold < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %g must be great than or equal zero",(double)jac->strongthreshold);
432:     PetscStackCallStandard(HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
433:   }
434:   PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);
435:   if (flg) {
436:     if (jac->maxrowsum < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %g must be greater than zero",(double)jac->maxrowsum);
437:     if (jac->maxrowsum > 1.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %g must be less than or equal one",(double)jac->maxrowsum);
438:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
439:   }

441:   /* Grid sweeps */
442:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_all","Number of sweeps for the up and down grid levels","None",jac->gridsweeps[0],&indx,&flg);
443:   if (flg) {
444:     PetscStackCallStandard(HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver,indx));
445:     /* modify the jac structure so we can view the updated options with PC_View */
446:     jac->gridsweeps[0] = indx;
447:     jac->gridsweeps[1] = indx;
448:     /*defaults coarse to 1 */
449:     jac->gridsweeps[2] = 1;
450:   }

452:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx,&flg);
453:   if (flg) {
454:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 1));
455:     jac->gridsweeps[0] = indx;
456:   }
457:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);
458:   if (flg) {
459:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 2));
460:     jac->gridsweeps[1] = indx;
461:   }
462:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);
463:   if (flg) {
464:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 3));
465:     jac->gridsweeps[2] = indx;
466:   }

468:   /* Relax type */
469:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_all","Relax type for the up and down cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);
470:   if (flg) {
471:     jac->relaxtype[0] = jac->relaxtype[1]  = indx;
472:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, indx));
473:     /* by default, coarse type set to 9 */
474:     jac->relaxtype[2] = 9;

476:   }
477:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);
478:   if (flg) {
479:     jac->relaxtype[0] = indx;
480:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 1));
481:   }
482:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);
483:   if (flg) {
484:     jac->relaxtype[1] = indx;
485:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 2));
486:   }
487:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[9],&indx,&flg);
488:   if (flg) {
489:     jac->relaxtype[2] = indx;
490:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 3));
491:   }

493:   /* Relaxation Weight */
494:   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);
495:   if (flg) {
496:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxWt,(jac->hsolver,tmpdbl));
497:     jac->relaxweight = tmpdbl;
498:   }

500:   n         = 2;
501:   twodbl[0] = twodbl[1] = 1.0;
502:   PetscOptionsRealArray("-pc_hypre_boomeramg_relax_weight_level","Set the relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);
503:   if (flg) {
504:     if (n == 2) {
505:       indx =  (int)PetscAbsReal(twodbl[1]);
506:       PetscStackCallStandard(HYPRE_BoomerAMGSetLevelRelaxWt,(jac->hsolver,twodbl[0],indx));
507:     } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relax weight level: you must provide 2 values separated by a comma (and no space), you provided %d",n);
508:   }

510:   /* Outer relaxation Weight */
511:   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);
512:   if (flg) {
513:     PetscStackCallStandard(HYPRE_BoomerAMGSetOuterWt,(jac->hsolver, tmpdbl));
514:     jac->outerrelaxweight = tmpdbl;
515:   }

517:   n         = 2;
518:   twodbl[0] = twodbl[1] = 1.0;
519:   PetscOptionsRealArray("-pc_hypre_boomeramg_outer_relax_weight_level","Set the outer relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);
520:   if (flg) {
521:     if (n == 2) {
522:       indx =  (int)PetscAbsReal(twodbl[1]);
523:       PetscStackCallStandard(HYPRE_BoomerAMGSetLevelOuterWt,(jac->hsolver, twodbl[0], indx));
524:     } else SETERRQ1(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 %d",n);
525:   }

527:   /* the Relax Order */
528:   PetscOptionsBool("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);

530:   if (flg && tmp_truth) {
531:     jac->relaxorder = 0;
532:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
533:   }
534:   PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,ALEN(HYPREBoomerAMGMeasureType),HYPREBoomerAMGMeasureType[0],&indx,&flg);
535:   if (flg) {
536:     jac->measuretype = indx;
537:     PetscStackCallStandard(HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
538:   }
539:   /* update list length 3/07 */
540:   PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,ALEN(HYPREBoomerAMGCoarsenType),HYPREBoomerAMGCoarsenType[6],&indx,&flg);
541:   if (flg) {
542:     jac->coarsentype = indx;
543:     PetscStackCallStandard(HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
544:   }

546:   /* new 3/07 */
547:   PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,ALEN(HYPREBoomerAMGInterpType),HYPREBoomerAMGInterpType[0],&indx,&flg);
548:   if (flg) {
549:     jac->interptype = indx;
550:     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
551:   }

553:   PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);
554:   if (flg) {
555:     level = 3;
556:     PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,NULL);

558:     jac->printstatistics = PETSC_TRUE;
559:     PetscStackCallStandard(HYPRE_BoomerAMGSetPrintLevel,(jac->hsolver,level));
560:   }

562:   PetscOptionsName("-pc_hypre_boomeramg_print_debug","Print debug information","None",&flg);
563:   if (flg) {
564:     level = 3;
565:     PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,NULL);

567:     jac->printstatistics = PETSC_TRUE;
568:     PetscStackCallStandard(HYPRE_BoomerAMGSetDebugFlag,(jac->hsolver,level));
569:   }

571:   PetscOptionsBool("-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", jac->nodal_coarsen ? PETSC_TRUE : PETSC_FALSE, &tmp_truth, &flg);
572:   if (flg && tmp_truth) {
573:     jac->nodal_coarsen = 1;
574:     PetscStackCallStandard(HYPRE_BoomerAMGSetNodal,(jac->hsolver,1));
575:   }

577:   PetscOptionsBool("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);
578:   if (flg && tmp_truth) {
579:     PetscInt tmp_int;
580:     PetscOptionsInt("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);
581:     if (flg) jac->nodal_relax_levels = tmp_int;
582:     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothType,(jac->hsolver,6));
583:     PetscStackCallStandard(HYPRE_BoomerAMGSetDomainType,(jac->hsolver,1));
584:     PetscStackCallStandard(HYPRE_BoomerAMGSetOverlap,(jac->hsolver,0));
585:     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothNumLevels,(jac->hsolver,jac->nodal_relax_levels));
586:   }

588:   PetscOptionsTail();
589:   return(0);
590: }

594: 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)
595: {
596:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
598:   PetscInt       oits;

601:   PetscCitationsRegister(hypreCitation,&cite);
602:   PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,its*jac->maxiter));
603:   PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,rtol));
604:   jac->applyrichardson = PETSC_TRUE;
605:   PCApply_HYPRE(pc,b,y);
606:   jac->applyrichardson = PETSC_FALSE;
607:   PetscStackCallStandard(HYPRE_BoomerAMGGetNumIterations,(jac->hsolver,(HYPRE_Int *)&oits));
608:   *outits = oits;
609:   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
610:   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
611:   PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
612:   PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
613:   return(0);
614: }


619: static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
620: {
621:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
623:   PetscBool      iascii;

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

639:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %g\n",(double)jac->maxrowsum);

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

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

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

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

670: /* --------------------------------------------------------------------------------------------*/
673: static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PetscOptions *PetscOptionsObject,PC pc)
674: {
675:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
677:   PetscInt       indx;
678:   PetscBool      flag;
679:   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};

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

687:   PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);
688:   if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter));

690:   PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);
691:   if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal));

693:   PetscOptionsBool("-pc_hypre_parasails_logging","Print info to screen","None",(PetscBool)jac->logging,(PetscBool*)&jac->logging,&flag);
694:   if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging));

696:   PetscOptionsBool("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscBool)jac->ruse,(PetscBool*)&jac->ruse,&flag);
697:   if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse));

699:   PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,ALEN(symtlist),symtlist[0],&indx,&flag);
700:   if (flag) {
701:     jac->symt = indx;
702:     PetscStackCallStandard(HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
703:   }

705:   PetscOptionsTail();
706:   return(0);
707: }

711: static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
712: {
713:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
715:   PetscBool      iascii;
716:   const char     *symt = 0;;

719:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
720:   if (iascii) {
721:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");
722:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);
723:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %g\n",(double)jac->threshhold);
724:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %g\n",(double)jac->filter);
725:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %g\n",(double)jac->loadbal);
726:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscBools[jac->ruse]);
727:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscBools[jac->logging]);
728:     if (!jac->symt) symt = "nonsymmetric matrix and preconditioner";
729:     else if (jac->symt == 1) symt = "SPD matrix and preconditioner";
730:     else if (jac->symt == 2) symt = "nonsymmetric matrix but SPD preconditioner";
731:     else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
732:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);
733:   }
734:   return(0);
735: }
736: /* --------------------------------------------------------------------------------------------*/
739: static PetscErrorCode PCSetFromOptions_HYPRE_AMS(PetscOptions *PetscOptionsObject,PC pc)
740: {
741:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
743:   PetscInt       n;
744:   PetscBool      flag,flag2,flag3,flag4;

747:   PetscOptionsHead(PetscOptionsObject,"HYPRE AMS Options");
748:   PetscOptionsInt("-pc_hypre_ams_print_level","Debugging output level for AMS","None",jac->as_print,&jac->as_print,&flag);
749:   if (flag) PetscStackCallStandard(HYPRE_AMSSetPrintLevel,(jac->hsolver,jac->as_print));
750:   PetscOptionsInt("-pc_hypre_ams_max_iter","Maximum number of AMS multigrid iterations within PCApply","None",jac->as_max_iter,&jac->as_max_iter,&flag);
751:   if (flag) PetscStackCallStandard(HYPRE_AMSSetMaxIter,(jac->hsolver,jac->as_max_iter));
752:   PetscOptionsInt("-pc_hypre_ams_cycle_type","Cycle type for AMS multigrid","None",jac->ams_cycle_type,&jac->ams_cycle_type,&flag);
753:   if (flag) PetscStackCallStandard(HYPRE_AMSSetCycleType,(jac->hsolver,jac->ams_cycle_type));
754:   PetscOptionsReal("-pc_hypre_ams_tol","Error tolerance for AMS multigrid","None",jac->as_tol,&jac->as_tol,&flag);
755:   if (flag) PetscStackCallStandard(HYPRE_AMSSetTol,(jac->hsolver,jac->as_tol));
756:   PetscOptionsInt("-pc_hypre_ams_relax_type","Relaxation type for AMS smoother","None",jac->as_relax_type,&jac->as_relax_type,&flag);
757:   PetscOptionsInt("-pc_hypre_ams_relax_times","Number of relaxation steps for AMS smoother","None",jac->as_relax_times,&jac->as_relax_times,&flag2);
758:   PetscOptionsReal("-pc_hypre_ams_relax_weight","Relaxation weight for AMS smoother","None",jac->as_relax_weight,&jac->as_relax_weight,&flag3);
759:   PetscOptionsReal("-pc_hypre_ams_omega","SSOR coefficient for AMS smoother","None",jac->as_omega,&jac->as_omega,&flag4);
760:   if (flag || flag2 || flag3 || flag4) {
761:     PetscStackCallStandard(HYPRE_AMSSetSmoothingOptions,(jac->hsolver,jac->as_relax_type,
762:                                                                       jac->as_relax_times,
763:                                                                       jac->as_relax_weight,
764:                                                                       jac->as_omega));
765:   }
766:   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);
767:   n = 5;
768:   PetscOptionsIntArray("-pc_hypre_ams_amg_alpha_options","AMG options for vector Poisson","None",jac->as_amg_alpha_opts,&n,&flag2);
769:   if (flag || flag2) {
770:     PetscStackCallStandard(HYPRE_AMSSetAlphaAMGOptions,(jac->hsolver,jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
771:                                                                      jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
772:                                                                      jac->as_amg_alpha_opts[2],       /* AMG relax_type */
773:                                                                      jac->as_amg_alpha_theta,
774:                                                                      jac->as_amg_alpha_opts[3],       /* AMG interp_type */
775:                                                                      jac->as_amg_alpha_opts[4]));     /* AMG Pmax */
776:   }
777:   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);
778:   n = 5;
779:   PetscOptionsIntArray("-pc_hypre_ams_amg_beta_options","AMG options for scalar Poisson solver","None",jac->as_amg_beta_opts,&n,&flag2);
780:   if (flag || flag2) {
781:     PetscStackCallStandard(HYPRE_AMSSetBetaAMGOptions,(jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
782:                                                                     jac->as_amg_beta_opts[1],       /* AMG agg_levels */
783:                                                                     jac->as_amg_beta_opts[2],       /* AMG relax_type */
784:                                                                     jac->as_amg_beta_theta,
785:                                                                     jac->as_amg_beta_opts[3],       /* AMG interp_type */
786:                                                                     jac->as_amg_beta_opts[4]));     /* AMG Pmax */
787:   }
788:   PetscOptionsTail();
789:   return(0);
790: }

794: static PetscErrorCode PCView_HYPRE_AMS(PC pc,PetscViewer viewer)
795: {
796:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
798:   PetscBool      iascii;

801:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
802:   if (iascii) {
803:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS preconditioning\n");
804:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: subspace iterations per application %d\n",jac->as_max_iter);
805:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: subspace cycle type %d\n",jac->ams_cycle_type);
806:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: subspace iteration tolerance %g\n",jac->as_tol);
807:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: smoother type %d\n",jac->as_relax_type);
808:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: number of smoothing steps %d\n",jac->as_relax_times);
809:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: smoother weight %g\n",jac->as_relax_weight);
810:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: smoother omega %g\n",jac->as_omega);
811:     if (jac->alpha_Poisson) {
812:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: vector Poisson solver (passed in by user)\n");
813:     } else {
814:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: vector Poisson solver (computed) \n");
815:     }
816:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG coarsening type %d\n",jac->as_amg_alpha_opts[0]);
817:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG levels of aggressive coarsening %d\n",jac->as_amg_alpha_opts[1]);
818:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG relaxation type %d\n",jac->as_amg_alpha_opts[2]);
819:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG interpolation type %d\n",jac->as_amg_alpha_opts[3]);
820:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG max nonzero elements in interpolation rows %d\n",jac->as_amg_alpha_opts[4]);
821:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG strength threshold %g\n",jac->as_amg_alpha_theta);
822:     if (!jac->ams_beta_is_zero) {
823:       if (jac->beta_Poisson) {
824:         PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: scalar Poisson solver (passed in by user)\n");
825:       } else {
826:         PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: scalar Poisson solver (computed) \n");
827:       }
828:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG coarsening type %d\n",jac->as_amg_beta_opts[0]);
829:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG levels of aggressive coarsening %d\n",jac->as_amg_beta_opts[1]);
830:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG relaxation type %d\n",jac->as_amg_beta_opts[2]);
831:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG interpolation type %d\n",jac->as_amg_beta_opts[3]);
832:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG max nonzero elements in interpolation rows %d\n",jac->as_amg_beta_opts[4]);
833:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG strength threshold %g\n",jac->as_amg_beta_theta);
834:     }
835:   }
836:   return(0);
837: }

841: static PetscErrorCode PCSetFromOptions_HYPRE_ADS(PetscOptions *PetscOptionsObject,PC pc)
842: {
843:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
845:   PetscInt       n;
846:   PetscBool      flag,flag2,flag3,flag4;

849:   PetscOptionsHead(PetscOptionsObject,"HYPRE ADS Options");
850:   PetscOptionsInt("-pc_hypre_ads_print_level","Debugging output level for ADS","None",jac->as_print,&jac->as_print,&flag);
851:   if (flag) PetscStackCallStandard(HYPRE_ADSSetPrintLevel,(jac->hsolver,jac->as_print));
852:   PetscOptionsInt("-pc_hypre_ads_max_iter","Maximum number of ADS multigrid iterations within PCApply","None",jac->as_max_iter,&jac->as_max_iter,&flag);
853:   if (flag) PetscStackCallStandard(HYPRE_ADSSetMaxIter,(jac->hsolver,jac->as_max_iter));
854:   PetscOptionsInt("-pc_hypre_ads_cycle_type","Cycle type for ADS multigrid","None",jac->ads_cycle_type,&jac->ads_cycle_type,&flag);
855:   if (flag) PetscStackCallStandard(HYPRE_ADSSetCycleType,(jac->hsolver,jac->ads_cycle_type));
856:   PetscOptionsReal("-pc_hypre_ads_tol","Error tolerance for ADS multigrid","None",jac->as_tol,&jac->as_tol,&flag);
857:   if (flag) PetscStackCallStandard(HYPRE_ADSSetTol,(jac->hsolver,jac->as_tol));
858:   PetscOptionsInt("-pc_hypre_ads_relax_type","Relaxation type for ADS smoother","None",jac->as_relax_type,&jac->as_relax_type,&flag);
859:   PetscOptionsInt("-pc_hypre_ads_relax_times","Number of relaxation steps for ADS smoother","None",jac->as_relax_times,&jac->as_relax_times,&flag2);
860:   PetscOptionsReal("-pc_hypre_ads_relax_weight","Relaxation weight for ADS smoother","None",jac->as_relax_weight,&jac->as_relax_weight,&flag3);
861:   PetscOptionsReal("-pc_hypre_ads_omega","SSOR coefficient for ADS smoother","None",jac->as_omega,&jac->as_omega,&flag4);
862:   if (flag || flag2 || flag3 || flag4) {
863:     PetscStackCallStandard(HYPRE_ADSSetSmoothingOptions,(jac->hsolver,jac->as_relax_type,
864:                                                                       jac->as_relax_times,
865:                                                                       jac->as_relax_weight,
866:                                                                       jac->as_omega));
867:   }
868:   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);
869:   n = 5;
870:   PetscOptionsIntArray("-pc_hypre_ads_ams_options","AMG options for AMS solver inside ADS","None",jac->as_amg_alpha_opts,&n,&flag2);
871:   PetscOptionsInt("-pc_hypre_ads_ams_cycle_type","Cycle type for AMS solver inside ADS","None",jac->ams_cycle_type,&jac->ams_cycle_type,&flag3);
872:   if (flag || flag2 || flag3) {
873:     PetscStackCallStandard(HYPRE_ADSSetAMSOptions,(jac->hsolver,jac->ams_cycle_type,             /* AMS cycle type */
874:                                                                 jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
875:                                                                 jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
876:                                                                 jac->as_amg_alpha_opts[2],       /* AMG relax_type */
877:                                                                 jac->as_amg_alpha_theta,
878:                                                                 jac->as_amg_alpha_opts[3],       /* AMG interp_type */
879:                                                                 jac->as_amg_alpha_opts[4]));     /* AMG Pmax */
880:   }
881:   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);
882:   n = 5;
883:   PetscOptionsIntArray("-pc_hypre_ads_amg_options","AMG options for vector AMG solver inside ADS","None",jac->as_amg_beta_opts,&n,&flag2);
884:   if (flag || flag2) {
885:     PetscStackCallStandard(HYPRE_ADSSetAMGOptions,(jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
886:                                                                 jac->as_amg_beta_opts[1],       /* AMG agg_levels */
887:                                                                 jac->as_amg_beta_opts[2],       /* AMG relax_type */
888:                                                                 jac->as_amg_beta_theta,
889:                                                                 jac->as_amg_beta_opts[3],       /* AMG interp_type */
890:                                                                 jac->as_amg_beta_opts[4]));     /* AMG Pmax */
891:   }
892:   PetscOptionsTail();
893:   return(0);
894: }

898: static PetscErrorCode PCView_HYPRE_ADS(PC pc,PetscViewer viewer)
899: {
900:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
902:   PetscBool      iascii;

905:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
906:   if (iascii) {
907:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS preconditioning\n");
908:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: subspace iterations per application %d\n",jac->as_max_iter);
909:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: subspace cycle type %d\n",jac->ads_cycle_type);
910:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: subspace iteration tolerance %g\n",jac->as_tol);
911:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: smoother type %d\n",jac->as_relax_type);
912:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: number of smoothing steps %d\n",jac->as_relax_times);
913:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: smoother weight %g\n",jac->as_relax_weight);
914:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: smoother omega %g\n",jac->as_omega);
915:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: AMS solver\n");
916:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     subspace cycle type %d\n",jac->ams_cycle_type);
917:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG coarsening type %d\n",jac->as_amg_alpha_opts[0]);
918:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG levels of aggressive coarsening %d\n",jac->as_amg_alpha_opts[1]);
919:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG relaxation type %d\n",jac->as_amg_alpha_opts[2]);
920:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG interpolation type %d\n",jac->as_amg_alpha_opts[3]);
921:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG max nonzero elements in interpolation rows %d\n",jac->as_amg_alpha_opts[4]);
922:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG strength threshold %g\n",jac->as_amg_alpha_theta);
923:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: vector Poisson solver\n");
924:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG coarsening type %d\n",jac->as_amg_beta_opts[0]);
925:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG levels of aggressive coarsening %d\n",jac->as_amg_beta_opts[1]);
926:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG relaxation type %d\n",jac->as_amg_beta_opts[2]);
927:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG interpolation type %d\n",jac->as_amg_beta_opts[3]);
928:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG max nonzero elements in interpolation rows %d\n",jac->as_amg_beta_opts[4]);
929:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG strength threshold %g\n",jac->as_amg_beta_theta);
930:   }
931:   return(0);
932: }

936: static PetscErrorCode PCHYPRESetDiscreteGradient_HYPRE(PC pc, Mat G)
937: {
938:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
939:   HYPRE_ParCSRMatrix parcsr_G;
940:   PetscErrorCode     ierr;

943:   /* throw away any discrete gradient if already set */
944:   if (jac->G) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->G));
945:   MatHYPRE_IJMatrixCreate(G,&jac->G);
946:   MatHYPRE_IJMatrixCopy(G,jac->G);
947:   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->G,(void**)(&parcsr_G)));
948:   PetscStackCall("Hypre set gradient",(*jac->setdgrad)(jac->hsolver,parcsr_G););
949:   return(0);
950: }

954: /*@
955:  PCHYPRESetDiscreteGradient - Set discrete gradient matrix

957:    Collective on PC

959:    Input Parameters:
960: +  pc - the preconditioning context
961: -  G - the discrete gradient

963:    Level: intermediate

965:    Notes: G should have as many rows as the number of edges and as many columns as the number of vertices in the mesh
966:           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

968: .seealso:
969: @*/
970: PetscErrorCode PCHYPRESetDiscreteGradient(PC pc, Mat G)
971: {

978:   PetscTryMethod(pc,"PCHYPRESetDiscreteGradient_C",(PC,Mat),(pc,G));
979:   return(0);
980: }

984: static PetscErrorCode PCHYPRESetDiscreteCurl_HYPRE(PC pc, Mat C)
985: {
986:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
987:   HYPRE_ParCSRMatrix parcsr_C;
988:   PetscErrorCode     ierr;

991:   /* throw away any discrete curl if already set */
992:   if (jac->C) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->C));
993:   MatHYPRE_IJMatrixCreate(C,&jac->C);
994:   MatHYPRE_IJMatrixCopy(C,jac->C);
995:   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->C,(void**)(&parcsr_C)));
996:   PetscStackCall("Hypre set curl",(*jac->setdcurl)(jac->hsolver,parcsr_C););
997:   return(0);
998: }

1002: /*@
1003:  PCHYPRESetDiscreteCurl - Set discrete curl matrix

1005:    Collective on PC

1007:    Input Parameters:
1008: +  pc - the preconditioning context
1009: -  C - the discrete curl

1011:    Level: intermediate

1013:    Notes: C should have as many rows as the number of faces and as many columns as the number of edges in the mesh
1014:           Each row of G has as many nonzeros as the number of edges of a face, with column indexes being the global indexes of the corresponding edge: matrix entries are +1 and -1 depending on edge orientation with respect to the face orientation

1016: .seealso:
1017: @*/
1018: PetscErrorCode PCHYPRESetDiscreteCurl(PC pc, Mat C)
1019: {

1026:   PetscTryMethod(pc,"PCHYPRESetDiscreteCurl_C",(PC,Mat),(pc,C));
1027:   return(0);
1028: }

1032: static PetscErrorCode PCHYPRESetAlphaPoissonMatrix_HYPRE_AMS(PC pc, Mat A)
1033: {
1034:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
1035:   HYPRE_ParCSRMatrix parcsr_alpha_Poisson;
1036:   PetscErrorCode     ierr;

1039:   /* throw away any matrix if already set */
1040:   if (jac->alpha_Poisson) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->alpha_Poisson));
1041:   MatHYPRE_IJMatrixCreate(A,&jac->alpha_Poisson);
1042:   MatHYPRE_IJMatrixCopy(A,jac->alpha_Poisson);
1043:   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->alpha_Poisson,(void**)(&parcsr_alpha_Poisson)));
1044:   PetscStackCallStandard(HYPRE_AMSSetAlphaPoissonMatrix,(jac->hsolver,parcsr_alpha_Poisson));
1045:   return(0);
1046: }

1050: /*@
1051:  PCHYPRESetAlphaPoissonMatrix - Set vector Poisson matrix

1053:    Collective on PC

1055:    Input Parameters:
1056: +  pc - the preconditioning context
1057: -  A - the matrix

1059:    Level: intermediate

1061:    Notes: A should be obtained by discretizing the vector valued Poisson problem with linear finite elements

1063: .seealso:
1064: @*/
1065: PetscErrorCode PCHYPRESetAlphaPoissonMatrix(PC pc, Mat A)
1066: {

1073:   PetscTryMethod(pc,"PCHYPRESetAlphaPoissonMatrix_C",(PC,Mat),(pc,A));
1074:   return(0);
1075: }

1079: static PetscErrorCode PCHYPRESetBetaPoissonMatrix_HYPRE_AMS(PC pc, Mat A)
1080: {
1081:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
1082:   HYPRE_ParCSRMatrix parcsr_beta_Poisson;
1083:   PetscErrorCode     ierr;

1086:   if (!A) {
1087:     PetscStackCallStandard(HYPRE_AMSSetBetaPoissonMatrix,(jac->hsolver,NULL));
1088:     jac->ams_beta_is_zero = PETSC_TRUE;
1089:     return(0);
1090:   }
1091:   jac->ams_beta_is_zero = PETSC_FALSE;
1092:   /* throw away any matrix if already set */
1093:   if (jac->beta_Poisson) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->beta_Poisson));
1094:   MatHYPRE_IJMatrixCreate(A,&jac->beta_Poisson);
1095:   MatHYPRE_IJMatrixCopy(A,jac->beta_Poisson);
1096:   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->beta_Poisson,(void**)(&parcsr_beta_Poisson)));
1097:   PetscStackCallStandard(HYPRE_AMSSetBetaPoissonMatrix,(jac->hsolver,parcsr_beta_Poisson));
1098:   return(0);
1099: }

1103: /*@
1104:  PCHYPRESetBetaPoissonMatrix - Set Poisson matrix

1106:    Collective on PC

1108:    Input Parameters:
1109: +  pc - the preconditioning context
1110: -  A - the matrix

1112:    Level: intermediate

1114:    Notes: A should be obtained by discretizing the Poisson problem with linear finite elements.
1115:           Following HYPRE convention, the scalar Poisson solver of AMS can be turned off by passing NULL.

1117: .seealso:
1118: @*/
1119: PetscErrorCode PCHYPRESetBetaPoissonMatrix(PC pc, Mat A)
1120: {

1125:   if (A) {
1128:   }
1129:   PetscTryMethod(pc,"PCHYPRESetBetaPoissonMatrix_C",(PC,Mat),(pc,A));
1130:   return(0);
1131: }

1135: static PetscErrorCode PCHYPRESetEdgeConstantVectors_HYPRE_AMS(PC pc,Vec ozz, Vec zoz, Vec zzo)
1136: {
1137:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
1138:   HYPRE_ParVector    par_ozz,par_zoz,par_zzo;
1139:   PetscInt           dim;
1140:   PetscErrorCode     ierr;

1143:   /* throw away any vector if already set */
1144:   if (jac->constants[0]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[0]));
1145:   if (jac->constants[1]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[1]));
1146:   if (jac->constants[2]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[2]));
1147:   jac->constants[0] = NULL;
1148:   jac->constants[1] = NULL;
1149:   jac->constants[2] = NULL;
1150:   VecHYPRE_IJVectorCreate(ozz,&jac->constants[0]);
1151:   VecHYPRE_IJVectorCopy(ozz,jac->constants[0]);
1152:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->constants[0],(void**)(&par_ozz)));
1153:   VecHYPRE_IJVectorCreate(zoz,&jac->constants[1]);
1154:   VecHYPRE_IJVectorCopy(zoz,jac->constants[1]);
1155:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->constants[1],(void**)(&par_zoz)));
1156:   dim = 2;
1157:   if (zzo) {
1158:     VecHYPRE_IJVectorCreate(zzo,&jac->constants[2]);
1159:     VecHYPRE_IJVectorCopy(zzo,jac->constants[2]);
1160:     PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->constants[2],(void**)(&par_zzo)));
1161:     dim++;
1162:   }
1163:   PetscStackCallStandard(HYPRE_AMSSetEdgeConstantVectors,(jac->hsolver,par_ozz,par_zoz,par_zzo));
1164:   PetscStackCallStandard(HYPRE_AMSSetDimension,(jac->hsolver,dim));
1165:   return(0);
1166: }

1170: /*@
1171:  PCHYPRESetEdgeConstantVectors - Set the representation of the constant vector fields in edge element basis

1173:    Collective on PC

1175:    Input Parameters:
1176: +  pc - the preconditioning context
1177: -  ozz - vector representing (1,0,0) (or (1,0) in 2D)
1178: -  zoz - vector representing (0,1,0) (or (0,1) in 2D)
1179: -  zzo - vector representing (0,0,1) (use NULL in 2D)

1181:    Level: intermediate

1183:    Notes:

1185: .seealso:
1186: @*/
1187: PetscErrorCode PCHYPRESetEdgeConstantVectors(PC pc, Vec ozz, Vec zoz, Vec zzo)
1188: {

1199:   PetscTryMethod(pc,"PCHYPRESetEdgeConstantVectors_C",(PC,Vec,Vec,Vec),(pc,ozz,zoz,zzo));
1200:   return(0);
1201: }

1205: static PetscErrorCode PCSetCoordinates_HYPRE(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords)
1206: {
1207:   PC_HYPRE        *jac = (PC_HYPRE*)pc->data;
1208:   Vec             tv;
1209:   HYPRE_ParVector par_coords[3];
1210:   PetscInt        i;
1211:   PetscErrorCode  ierr;

1214:   /* throw away any coordinate vector if already set */
1215:   if (jac->coords[0]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[0]));
1216:   if (jac->coords[1]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[1]));
1217:   if (jac->coords[2]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[2]));
1218:   /* set problem's dimension */
1219:   if (jac->setdim) {
1220:     PetscStackCall("Hypre set dim",(*jac->setdim)(jac->hsolver,dim););
1221:   }
1222:   /* compute IJ vector for coordinates */
1223:   VecCreate(PetscObjectComm((PetscObject)pc),&tv);
1224:   VecSetType(tv,VECSTANDARD);
1225:   VecSetSizes(tv,nloc,PETSC_DECIDE);
1226:   for (i=0;i<dim;i++) {
1227:     PetscScalar *array;
1228:     PetscInt    j;

1230:     VecHYPRE_IJVectorCreate(tv,&jac->coords[i]);
1231:     VecGetArray(tv,&array);
1232:     for (j=0;j<nloc;j++) {
1233:       array[j] = coords[j*dim+i];
1234:     }
1235:     PetscStackCallStandard(HYPRE_IJVectorSetValues,(jac->coords[i],nloc,NULL,array));
1236:     PetscStackCallStandard(HYPRE_IJVectorAssemble,(jac->coords[i]));
1237:     VecRestoreArray(tv,&array);
1238:   }
1239:   VecDestroy(&tv);
1240:   /* pass parCSR vectors to AMS solver */
1241:   par_coords[0] = NULL;
1242:   par_coords[1] = NULL;
1243:   par_coords[2] = NULL;
1244:   if (jac->coords[0]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[0],(void**)(&par_coords[0])));
1245:   if (jac->coords[1]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[1],(void**)(&par_coords[1])));
1246:   if (jac->coords[2]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[2],(void**)(&par_coords[2])));
1247:   PetscStackCall("Hypre set coords",(*jac->setcoord)(jac->hsolver,par_coords[0],par_coords[1],par_coords[2]););
1248:   return(0);
1249: }

1251: /* ---------------------------------------------------------------------------------*/

1255: static PetscErrorCode  PCHYPREGetType_HYPRE(PC pc,const char *name[])
1256: {
1257:   PC_HYPRE *jac = (PC_HYPRE*)pc->data;

1260:   *name = jac->hypre_type;
1261:   return(0);
1262: }

1266: static PetscErrorCode  PCHYPRESetType_HYPRE(PC pc,const char name[])
1267: {
1268:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1270:   PetscBool      flag;

1273:   if (jac->hypre_type) {
1274:     PetscStrcmp(jac->hypre_type,name,&flag);
1275:     if (!flag) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
1276:     return(0);
1277:   } else {
1278:     PetscStrallocpy(name, &jac->hypre_type);
1279:   }

1281:   jac->maxiter         = PETSC_DEFAULT;
1282:   jac->tol             = PETSC_DEFAULT;
1283:   jac->printstatistics = PetscLogPrintInfo;

1285:   PetscStrcmp("pilut",jac->hypre_type,&flag);
1286:   if (flag) {
1287:     PetscStackCallStandard(HYPRE_ParCSRPilutCreate,(jac->comm_hypre,&jac->hsolver));
1288:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
1289:     pc->ops->view           = PCView_HYPRE_Pilut;
1290:     jac->destroy            = HYPRE_ParCSRPilutDestroy;
1291:     jac->setup              = HYPRE_ParCSRPilutSetup;
1292:     jac->solve              = HYPRE_ParCSRPilutSolve;
1293:     jac->factorrowsize      = PETSC_DEFAULT;
1294:     return(0);
1295:   }
1296:   PetscStrcmp("parasails",jac->hypre_type,&flag);
1297:   if (flag) {
1298:     PetscStackCallStandard(HYPRE_ParaSailsCreate,(jac->comm_hypre,&jac->hsolver));
1299:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
1300:     pc->ops->view           = PCView_HYPRE_ParaSails;
1301:     jac->destroy            = HYPRE_ParaSailsDestroy;
1302:     jac->setup              = HYPRE_ParaSailsSetup;
1303:     jac->solve              = HYPRE_ParaSailsSolve;
1304:     /* initialize */
1305:     jac->nlevels    = 1;
1306:     jac->threshhold = .1;
1307:     jac->filter     = .1;
1308:     jac->loadbal    = 0;
1309:     if (PetscLogPrintInfo) jac->logging = (int) PETSC_TRUE;
1310:     else jac->logging = (int) PETSC_FALSE;

1312:     jac->ruse = (int) PETSC_FALSE;
1313:     jac->symt = 0;
1314:     PetscStackCallStandard(HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels));
1315:     PetscStackCallStandard(HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter));
1316:     PetscStackCallStandard(HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal));
1317:     PetscStackCallStandard(HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging));
1318:     PetscStackCallStandard(HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse));
1319:     PetscStackCallStandard(HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
1320:     return(0);
1321:   }
1322:   PetscStrcmp("boomeramg",jac->hypre_type,&flag);
1323:   if (flag) {
1324:     HYPRE_BoomerAMGCreate(&jac->hsolver);
1325:     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
1326:     pc->ops->view            = PCView_HYPRE_BoomerAMG;
1327:     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
1328:     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
1329:     jac->destroy             = HYPRE_BoomerAMGDestroy;
1330:     jac->setup               = HYPRE_BoomerAMGSetup;
1331:     jac->solve               = HYPRE_BoomerAMGSolve;
1332:     jac->applyrichardson     = PETSC_FALSE;
1333:     /* these defaults match the hypre defaults */
1334:     jac->cycletype        = 1;
1335:     jac->maxlevels        = 25;
1336:     jac->maxiter          = 1;
1337:     jac->tol              = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */
1338:     jac->truncfactor      = 0.0;
1339:     jac->strongthreshold  = .25;
1340:     jac->maxrowsum        = .9;
1341:     jac->coarsentype      = 6;
1342:     jac->measuretype      = 0;
1343:     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = 1;
1344:     jac->relaxtype[0]     = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */
1345:     jac->relaxtype[2]     = 9; /*G.E. */
1346:     jac->relaxweight      = 1.0;
1347:     jac->outerrelaxweight = 1.0;
1348:     jac->relaxorder       = 1;
1349:     jac->interptype       = 0;
1350:     jac->agg_nl           = 0;
1351:     jac->pmax             = 0;
1352:     jac->truncfactor      = 0.0;
1353:     jac->agg_num_paths    = 1;

1355:     jac->nodal_coarsen      = 0;
1356:     jac->nodal_relax        = PETSC_FALSE;
1357:     jac->nodal_relax_levels = 1;
1358:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
1359:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
1360:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
1361:     PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
1362:     PetscStackCallStandard(HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
1363:     PetscStackCallStandard(HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
1364:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
1365:     PetscStackCallStandard(HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
1366:     PetscStackCallStandard(HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
1367:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
1368:     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
1369:     PetscStackCallStandard(HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
1370:     PetscStackCallStandard(HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
1371:     PetscStackCallStandard(HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
1372:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, jac->relaxtype[0]));  /*defaults coarse to 9*/
1373:     PetscStackCallStandard(HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver, jac->gridsweeps[0])); /*defaults coarse to 1 */
1374:     return(0);
1375:   }
1376:   PetscStrcmp("ams",jac->hypre_type,&flag);
1377:   if (flag) {
1378:     HYPRE_AMSCreate(&jac->hsolver);
1379:     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_AMS;
1380:     pc->ops->view            = PCView_HYPRE_AMS;
1381:     jac->destroy             = HYPRE_AMSDestroy;
1382:     jac->setup               = HYPRE_AMSSetup;
1383:     jac->solve               = HYPRE_AMSSolve;
1384:     jac->setdgrad            = HYPRE_AMSSetDiscreteGradient;
1385:     jac->setcoord            = HYPRE_AMSSetCoordinateVectors;
1386:     jac->setdim              = HYPRE_AMSSetDimension;
1387:     jac->coords[0]           = NULL;
1388:     jac->coords[1]           = NULL;
1389:     jac->coords[2]           = NULL;
1390:     jac->G                   = NULL;
1391:     /* solver parameters: these are borrowed from mfem package, and they are not the default values from HYPRE AMS */
1392:     jac->as_print           = 0;
1393:     jac->as_max_iter        = 1; /* used as a preconditioner */
1394:     jac->as_tol             = 0.; /* used as a preconditioner */
1395:     jac->ams_cycle_type     = 13;
1396:     /* Smoothing options */
1397:     jac->as_relax_type      = 2;
1398:     jac->as_relax_times     = 1;
1399:     jac->as_relax_weight    = 1.0;
1400:     jac->as_omega           = 1.0;
1401:     /* Vector valued Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
1402:     jac->as_amg_alpha_opts[0] = 10;
1403:     jac->as_amg_alpha_opts[1] = 1;
1404:     jac->as_amg_alpha_opts[2] = 6;
1405:     jac->as_amg_alpha_opts[3] = 6;
1406:     jac->as_amg_alpha_opts[4] = 4;
1407:     jac->as_amg_alpha_theta   = 0.25;
1408:     /* Scalar Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
1409:     jac->ams_beta_is_zero = PETSC_FALSE;
1410:     jac->as_amg_beta_opts[0] = 10;
1411:     jac->as_amg_beta_opts[1] = 1;
1412:     jac->as_amg_beta_opts[2] = 6;
1413:     jac->as_amg_beta_opts[3] = 6;
1414:     jac->as_amg_beta_opts[4] = 4;
1415:     jac->as_amg_beta_theta   = 0.25;
1416:     PetscStackCallStandard(HYPRE_AMSSetPrintLevel,(jac->hsolver,jac->as_print));
1417:     PetscStackCallStandard(HYPRE_AMSSetMaxIter,(jac->hsolver,jac->as_max_iter));
1418:     PetscStackCallStandard(HYPRE_AMSSetCycleType,(jac->hsolver,jac->ams_cycle_type));
1419:     PetscStackCallStandard(HYPRE_AMSSetTol,(jac->hsolver,jac->as_tol));
1420:     PetscStackCallStandard(HYPRE_AMSSetSmoothingOptions,(jac->hsolver,jac->as_relax_type,
1421:                                                                       jac->as_relax_times,
1422:                                                                       jac->as_relax_weight,
1423:                                                                       jac->as_omega));
1424:     PetscStackCallStandard(HYPRE_AMSSetAlphaAMGOptions,(jac->hsolver,jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
1425:                                                                      jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
1426:                                                                      jac->as_amg_alpha_opts[2],       /* AMG relax_type */
1427:                                                                      jac->as_amg_alpha_theta,
1428:                                                                      jac->as_amg_alpha_opts[3],       /* AMG interp_type */
1429:                                                                      jac->as_amg_alpha_opts[4]));     /* AMG Pmax */
1430:     PetscStackCallStandard(HYPRE_AMSSetBetaAMGOptions,(jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
1431:                                                                     jac->as_amg_beta_opts[1],       /* AMG agg_levels */
1432:                                                                     jac->as_amg_beta_opts[2],       /* AMG relax_type */
1433:                                                                     jac->as_amg_beta_theta,
1434:                                                                     jac->as_amg_beta_opts[3],       /* AMG interp_type */
1435:                                                                     jac->as_amg_beta_opts[4]));     /* AMG Pmax */
1436:     PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_HYPRE);
1437:     PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteGradient_C",PCHYPRESetDiscreteGradient_HYPRE);
1438:     PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetEdgeConstantVectors_C",PCHYPRESetEdgeConstantVectors_HYPRE_AMS);
1439:     PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetAlphaPoissonMatrix_C",PCHYPRESetAlphaPoissonMatrix_HYPRE_AMS);
1440:     PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetBetaPoissonMatrix_C",PCHYPRESetBetaPoissonMatrix_HYPRE_AMS);
1441:     return(0);
1442:   }
1443:   PetscStrcmp("ads",jac->hypre_type,&flag);
1444:   if (flag) {
1445:     HYPRE_ADSCreate(&jac->hsolver);
1446:     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_ADS;
1447:     pc->ops->view            = PCView_HYPRE_ADS;
1448:     jac->destroy             = HYPRE_ADSDestroy;
1449:     jac->setup               = HYPRE_ADSSetup;
1450:     jac->solve               = HYPRE_ADSSolve;
1451:     jac->setdgrad            = HYPRE_ADSSetDiscreteGradient;
1452:     jac->setdcurl            = HYPRE_ADSSetDiscreteCurl;
1453:     jac->setcoord            = HYPRE_ADSSetCoordinateVectors;
1454:     jac->coords[0]           = NULL;
1455:     jac->coords[1]           = NULL;
1456:     jac->coords[2]           = NULL;
1457:     jac->G                   = NULL;
1458:     jac->C                   = NULL;
1459:     /* solver parameters: these are borrowed from mfem package, and they are not the default values from HYPRE ADS */
1460:     jac->as_print           = 0;
1461:     jac->as_max_iter        = 1; /* used as a preconditioner */
1462:     jac->as_tol             = 0.; /* used as a preconditioner */
1463:     jac->ads_cycle_type     = 13;
1464:     /* Smoothing options */
1465:     jac->as_relax_type      = 2;
1466:     jac->as_relax_times     = 1;
1467:     jac->as_relax_weight    = 1.0;
1468:     jac->as_omega           = 1.0;
1469:     /* AMS solver parameters: cycle_type, coarsen type, agg_levels, relax_type, interp_type, Pmax */
1470:     jac->ams_cycle_type       = 14;
1471:     jac->as_amg_alpha_opts[0] = 10;
1472:     jac->as_amg_alpha_opts[1] = 1;
1473:     jac->as_amg_alpha_opts[2] = 6;
1474:     jac->as_amg_alpha_opts[3] = 6;
1475:     jac->as_amg_alpha_opts[4] = 4;
1476:     jac->as_amg_alpha_theta   = 0.25;
1477:     /* Vector Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
1478:     jac->as_amg_beta_opts[0] = 10;
1479:     jac->as_amg_beta_opts[1] = 1;
1480:     jac->as_amg_beta_opts[2] = 6;
1481:     jac->as_amg_beta_opts[3] = 6;
1482:     jac->as_amg_beta_opts[4] = 4;
1483:     jac->as_amg_beta_theta   = 0.25;
1484:     PetscStackCallStandard(HYPRE_ADSSetPrintLevel,(jac->hsolver,jac->as_print));
1485:     PetscStackCallStandard(HYPRE_ADSSetMaxIter,(jac->hsolver,jac->as_max_iter));
1486:     PetscStackCallStandard(HYPRE_ADSSetCycleType,(jac->hsolver,jac->ams_cycle_type));
1487:     PetscStackCallStandard(HYPRE_ADSSetTol,(jac->hsolver,jac->as_tol));
1488:     PetscStackCallStandard(HYPRE_ADSSetSmoothingOptions,(jac->hsolver,jac->as_relax_type,
1489:                                                                       jac->as_relax_times,
1490:                                                                       jac->as_relax_weight,
1491:                                                                       jac->as_omega));
1492:     PetscStackCallStandard(HYPRE_ADSSetAMSOptions,(jac->hsolver,jac->ams_cycle_type,             /* AMG coarsen type */
1493:                                                                 jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
1494:                                                                 jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
1495:                                                                 jac->as_amg_alpha_opts[2],       /* AMG relax_type */
1496:                                                                 jac->as_amg_alpha_theta,
1497:                                                                 jac->as_amg_alpha_opts[3],       /* AMG interp_type */
1498:                                                                 jac->as_amg_alpha_opts[4]));     /* AMG Pmax */
1499:     PetscStackCallStandard(HYPRE_ADSSetAMGOptions,(jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
1500:                                                                 jac->as_amg_beta_opts[1],       /* AMG agg_levels */
1501:                                                                 jac->as_amg_beta_opts[2],       /* AMG relax_type */
1502:                                                                 jac->as_amg_beta_theta,
1503:                                                                 jac->as_amg_beta_opts[3],       /* AMG interp_type */
1504:                                                                 jac->as_amg_beta_opts[4]));     /* AMG Pmax */
1505:     PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_HYPRE);
1506:     PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteGradient_C",PCHYPRESetDiscreteGradient_HYPRE);
1507:     PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteCurl_C",PCHYPRESetDiscreteCurl_HYPRE);
1508:     return(0);
1509:   }
1510:   PetscFree(jac->hypre_type);

1512:   jac->hypre_type = NULL;
1513:   SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, boomeramg, ams",name);
1514:   return(0);
1515: }

1517: /*
1518:     It only gets here if the HYPRE type has not been set before the call to
1519:    ...SetFromOptions() which actually is most of the time
1520: */
1523: static PetscErrorCode PCSetFromOptions_HYPRE(PetscOptions *PetscOptionsObject,PC pc)
1524: {
1526:   PetscInt       indx;
1527:   const char     *type[] = {"pilut","parasails","boomeramg","ams","ads"};
1528:   PetscBool      flg;

1531:   PetscOptionsHead(PetscOptionsObject,"HYPRE preconditioner options");
1532:   PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"boomeramg",&indx,&flg);
1533:   if (flg) {
1534:     PCHYPRESetType_HYPRE(pc,type[indx]);
1535:   } else {
1536:     PCHYPRESetType_HYPRE(pc,"boomeramg");
1537:   }
1538:   if (pc->ops->setfromoptions) {
1539:     pc->ops->setfromoptions(PetscOptionsObject,pc);
1540:   }
1541:   PetscOptionsTail();
1542:   return(0);
1543: }

1547: /*@C
1548:      PCHYPRESetType - Sets which hypre preconditioner you wish to use

1550:    Input Parameters:
1551: +     pc - the preconditioner context
1552: -     name - either  pilut, parasails, boomeramg, ams, ads

1554:    Options Database Keys:
1555:    -pc_hypre_type - One of pilut, parasails, boomeramg, ams, ads

1557:    Level: intermediate

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

1562: @*/
1563: PetscErrorCode  PCHYPRESetType(PC pc,const char name[])
1564: {

1570:   PetscTryMethod(pc,"PCHYPRESetType_C",(PC,const char[]),(pc,name));
1571:   return(0);
1572: }

1576: /*@C
1577:      PCHYPREGetType - Gets which hypre preconditioner you are using

1579:    Input Parameter:
1580: .     pc - the preconditioner context

1582:    Output Parameter:
1583: .     name - either  pilut, parasails, boomeramg, ams, ads

1585:    Level: intermediate

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

1590: @*/
1591: PetscErrorCode  PCHYPREGetType(PC pc,const char *name[])
1592: {

1598:   PetscTryMethod(pc,"PCHYPREGetType_C",(PC,const char*[]),(pc,name));
1599:   return(0);
1600: }

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

1605:    Options Database Keys:
1606: +   -pc_hypre_type - One of pilut, parasails, boomeramg, ams, ads
1607: -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
1608:           preconditioner

1610:    Level: intermediate

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

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

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

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

1634:           See PCPFMG for access to the hypre Struct PFMG solver

1636: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
1637:            PCHYPRESetType(), PCPFMG

1639: M*/

1643: PETSC_EXTERN PetscErrorCode PCCreate_HYPRE(PC pc)
1644: {
1645:   PC_HYPRE       *jac;

1649:   PetscNewLog(pc,&jac);

1651:   pc->data                = jac;
1652:   pc->ops->destroy        = PCDestroy_HYPRE;
1653:   pc->ops->setfromoptions = PCSetFromOptions_HYPRE;
1654:   pc->ops->setup          = PCSetUp_HYPRE;
1655:   pc->ops->apply          = PCApply_HYPRE;
1656:   jac->comm_hypre         = MPI_COMM_NULL;
1657:   jac->hypre_type         = NULL;
1658:   jac->coords[0]          = NULL;
1659:   jac->coords[1]          = NULL;
1660:   jac->coords[2]          = NULL;
1661:   jac->constants[0]       = NULL;
1662:   jac->constants[1]       = NULL;
1663:   jac->constants[2]       = NULL;
1664:   jac->G                  = NULL;
1665:   jac->C                  = NULL;
1666:   jac->alpha_Poisson      = NULL;
1667:   jac->beta_Poisson       = NULL;
1668:   jac->setdim             = NULL;
1669:   /* duplicate communicator for hypre */
1670:   MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(jac->comm_hypre));
1671:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetType_C",PCHYPRESetType_HYPRE);
1672:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPREGetType_C",PCHYPREGetType_HYPRE);
1673:   return(0);
1674: }

1676: /* ---------------------------------------------------------------------------------------------------------------------------------*/

1678: /* this include is needed ONLY to allow access to the private data inside the Mat object specific to hypre */
1679: #include <petsc/private/matimpl.h>

1681: typedef struct {
1682:   MPI_Comm           hcomm;        /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */
1683:   HYPRE_StructSolver hsolver;

1685:   /* keep copy of PFMG options used so may view them */
1686:   PetscInt its;
1687:   double   tol;
1688:   PetscInt relax_type;
1689:   PetscInt rap_type;
1690:   PetscInt num_pre_relax,num_post_relax;
1691:   PetscInt max_levels;
1692: } PC_PFMG;

1696: PetscErrorCode PCDestroy_PFMG(PC pc)
1697: {
1699:   PC_PFMG        *ex = (PC_PFMG*) pc->data;

1702:   if (ex->hsolver) PetscStackCallStandard(HYPRE_StructPFMGDestroy,(ex->hsolver));
1703:   MPI_Comm_free(&ex->hcomm);
1704:   PetscFree(pc->data);
1705:   return(0);
1706: }

1708: static const char *PFMGRelaxType[] = {"Jacobi","Weighted-Jacobi","symmetric-Red/Black-Gauss-Seidel","Red/Black-Gauss-Seidel"};
1709: static const char *PFMGRAPType[] = {"Galerkin","non-Galerkin"};

1713: PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer)
1714: {
1716:   PetscBool      iascii;
1717:   PC_PFMG        *ex = (PC_PFMG*) pc->data;

1720:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1721:   if (iascii) {
1722:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG preconditioning\n");
1723:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max iterations %d\n",ex->its);
1724:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: tolerance %g\n",ex->tol);
1725:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);
1726:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: RAP type %s\n",PFMGRAPType[ex->rap_type]);
1727:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);
1728:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max levels %d\n",ex->max_levels);
1729:   }
1730:   return(0);
1731: }


1736: PetscErrorCode PCSetFromOptions_PFMG(PetscOptions *PetscOptionsObject,PC pc)
1737: {
1739:   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1740:   PetscBool      flg = PETSC_FALSE;

1743:   PetscOptionsHead(PetscOptionsObject,"PFMG options");
1744:   PetscOptionsBool("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",flg,&flg,NULL);
1745:   if (flg) {
1746:     int level=3;
1747:     PetscStackCallStandard(HYPRE_StructPFMGSetPrintLevel,(ex->hsolver,level));
1748:   }
1749:   PetscOptionsInt("-pc_pfmg_its","Number of iterations of PFMG to use as preconditioner","HYPRE_StructPFMGSetMaxIter",ex->its,&ex->its,NULL);
1750:   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(ex->hsolver,ex->its));
1751:   PetscOptionsInt("-pc_pfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_StructPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,NULL);
1752:   PetscStackCallStandard(HYPRE_StructPFMGSetNumPreRelax,(ex->hsolver,ex->num_pre_relax));
1753:   PetscOptionsInt("-pc_pfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_StructPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,NULL);
1754:   PetscStackCallStandard(HYPRE_StructPFMGSetNumPostRelax,(ex->hsolver,ex->num_post_relax));

1756:   PetscOptionsInt("-pc_pfmg_max_levels","Max Levels for MG hierarchy","HYPRE_StructPFMGSetMaxLevels",ex->max_levels,&ex->max_levels,NULL);
1757:   PetscStackCallStandard(HYPRE_StructPFMGSetMaxLevels,(ex->hsolver,ex->max_levels));

1759:   PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,NULL);
1760:   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(ex->hsolver,ex->tol));
1761:   PetscOptionsEList("-pc_pfmg_relax_type","Relax type for the up and down cycles","HYPRE_StructPFMGSetRelaxType",PFMGRelaxType,ALEN(PFMGRelaxType),PFMGRelaxType[ex->relax_type],&ex->relax_type,NULL);
1762:   PetscStackCallStandard(HYPRE_StructPFMGSetRelaxType,(ex->hsolver, ex->relax_type));
1763:   PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,ALEN(PFMGRAPType),PFMGRAPType[ex->rap_type],&ex->rap_type,NULL);
1764:   PetscStackCallStandard(HYPRE_StructPFMGSetRAPType,(ex->hsolver, ex->rap_type));
1765:   PetscOptionsTail();
1766:   return(0);
1767: }

1771: PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y)
1772: {
1773:   PetscErrorCode    ierr;
1774:   PC_PFMG           *ex = (PC_PFMG*) pc->data;
1775:   PetscScalar       *yy;
1776:   const PetscScalar *xx;
1777:   PetscInt          ilower[3],iupper[3];
1778:   Mat_HYPREStruct   *mx = (Mat_HYPREStruct*)(pc->pmat->data);

1781:   PetscCitationsRegister(hypreCitation,&cite);
1782:   DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);
1783:   iupper[0] += ilower[0] - 1;
1784:   iupper[1] += ilower[1] - 1;
1785:   iupper[2] += ilower[2] - 1;

1787:   /* copy x values over to hypre */
1788:   PetscStackCallStandard(HYPRE_StructVectorSetConstantValues,(mx->hb,0.0));
1789:   VecGetArrayRead(x,&xx);
1790:   PetscStackCallStandard(HYPRE_StructVectorSetBoxValues,(mx->hb,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,(PetscScalar*)xx));
1791:   VecRestoreArrayRead(x,&xx);
1792:   PetscStackCallStandard(HYPRE_StructVectorAssemble,(mx->hb));
1793:   PetscStackCallStandard(HYPRE_StructPFMGSolve,(ex->hsolver,mx->hmat,mx->hb,mx->hx));

1795:   /* copy solution values back to PETSc */
1796:   VecGetArray(y,&yy);
1797:   PetscStackCallStandard(HYPRE_StructVectorGetBoxValues,(mx->hx,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,yy));
1798:   VecRestoreArray(y,&yy);
1799:   return(0);
1800: }

1804: 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)
1805: {
1806:   PC_PFMG        *jac = (PC_PFMG*)pc->data;
1808:   PetscInt       oits;

1811:   PetscCitationsRegister(hypreCitation,&cite);
1812:   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(jac->hsolver,its*jac->its));
1813:   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(jac->hsolver,rtol));

1815:   PCApply_PFMG(pc,b,y);
1816:   PetscStackCallStandard(HYPRE_StructPFMGGetNumIterations,(jac->hsolver,(HYPRE_Int *)&oits));
1817:   *outits = oits;
1818:   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
1819:   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
1820:   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(jac->hsolver,jac->tol));
1821:   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(jac->hsolver,jac->its));
1822:   return(0);
1823: }


1828: PetscErrorCode PCSetUp_PFMG(PC pc)
1829: {
1830:   PetscErrorCode  ierr;
1831:   PC_PFMG         *ex = (PC_PFMG*) pc->data;
1832:   Mat_HYPREStruct *mx = (Mat_HYPREStruct*)(pc->pmat->data);
1833:   PetscBool       flg;

1836:   PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESTRUCT,&flg);
1837:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESTRUCT with this preconditioner");

1839:   /* create the hypre solver object and set its information */
1840:   if (ex->hsolver) PetscStackCallStandard(HYPRE_StructPFMGDestroy,(ex->hsolver));
1841:   PetscStackCallStandard(HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
1842:   PetscStackCallStandard(HYPRE_StructPFMGSetup,(ex->hsolver,mx->hmat,mx->hb,mx->hx));
1843:   PetscStackCallStandard(HYPRE_StructPFMGSetZeroGuess,(ex->hsolver));
1844:   return(0);
1845: }


1848: /*MC
1849:      PCPFMG - the hypre PFMG multigrid solver

1851:    Level: advanced

1853:    Options Database:
1854: + -pc_pfmg_its <its> number of iterations of PFMG to use as preconditioner
1855: . -pc_pfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
1856: . -pc_pfmg_num_post_relax <steps> number of smoothing steps after coarse grid
1857: . -pc_pfmg_tol <tol> tolerance of PFMG
1858: . -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
1859: - -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin

1861:    Notes:  This is for CELL-centered descretizations

1863:            This must be used with the MATHYPRESTRUCT matrix type.
1864:            This is less general than in hypre, it supports only one block per process defined by a PETSc DMDA.

1866: .seealso:  PCMG, MATHYPRESTRUCT
1867: M*/

1871: PETSC_EXTERN PetscErrorCode PCCreate_PFMG(PC pc)
1872: {
1874:   PC_PFMG        *ex;

1877:   PetscNew(&ex); \
1878:   pc->data = ex;

1880:   ex->its            = 1;
1881:   ex->tol            = 1.e-8;
1882:   ex->relax_type     = 1;
1883:   ex->rap_type       = 0;
1884:   ex->num_pre_relax  = 1;
1885:   ex->num_post_relax = 1;
1886:   ex->max_levels     = 0;

1888:   pc->ops->setfromoptions  = PCSetFromOptions_PFMG;
1889:   pc->ops->view            = PCView_PFMG;
1890:   pc->ops->destroy         = PCDestroy_PFMG;
1891:   pc->ops->apply           = PCApply_PFMG;
1892:   pc->ops->applyrichardson = PCApplyRichardson_PFMG;
1893:   pc->ops->setup           = PCSetUp_PFMG;

1895:   MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ex->hcomm));
1896:   PetscStackCallStandard(HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
1897:   return(0);
1898: }

1900: /* ---------------------------------------------------------------------------------------------------------------------------------------------------*/

1902: /* we know we are working with a HYPRE_SStructMatrix */
1903: typedef struct {
1904:   MPI_Comm            hcomm;       /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */
1905:   HYPRE_SStructSolver ss_solver;

1907:   /* keep copy of SYSPFMG options used so may view them */
1908:   PetscInt its;
1909:   double   tol;
1910:   PetscInt relax_type;
1911:   PetscInt num_pre_relax,num_post_relax;
1912: } PC_SysPFMG;

1916: PetscErrorCode PCDestroy_SysPFMG(PC pc)
1917: {
1919:   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;

1922:   if (ex->ss_solver) PetscStackCallStandard(HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));
1923:   MPI_Comm_free(&ex->hcomm);
1924:   PetscFree(pc->data);
1925:   return(0);
1926: }

1928: static const char *SysPFMGRelaxType[] = {"Weighted-Jacobi","Red/Black-Gauss-Seidel"};

1932: PetscErrorCode PCView_SysPFMG(PC pc,PetscViewer viewer)
1933: {
1935:   PetscBool      iascii;
1936:   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;

1939:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1940:   if (iascii) {
1941:     PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG preconditioning\n");
1942:     PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: max iterations %d\n",ex->its);
1943:     PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: tolerance %g\n",ex->tol);
1944:     PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);
1945:     PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);
1946:   }
1947:   return(0);
1948: }


1953: PetscErrorCode PCSetFromOptions_SysPFMG(PetscOptions *PetscOptionsObject,PC pc)
1954: {
1956:   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;
1957:   PetscBool      flg = PETSC_FALSE;

1960:   PetscOptionsHead(PetscOptionsObject,"SysPFMG options");
1961:   PetscOptionsBool("-pc_syspfmg_print_statistics","Print statistics","HYPRE_SStructSysPFMGSetPrintLevel",flg,&flg,NULL);
1962:   if (flg) {
1963:     int level=3;
1964:     PetscStackCallStandard(HYPRE_SStructSysPFMGSetPrintLevel,(ex->ss_solver,level));
1965:   }
1966:   PetscOptionsInt("-pc_syspfmg_its","Number of iterations of SysPFMG to use as preconditioner","HYPRE_SStructSysPFMGSetMaxIter",ex->its,&ex->its,NULL);
1967:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(ex->ss_solver,ex->its));
1968:   PetscOptionsInt("-pc_syspfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_SStructSysPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,NULL);
1969:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetNumPreRelax,(ex->ss_solver,ex->num_pre_relax));
1970:   PetscOptionsInt("-pc_syspfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_SStructSysPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,NULL);
1971:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetNumPostRelax,(ex->ss_solver,ex->num_post_relax));

1973:   PetscOptionsReal("-pc_syspfmg_tol","Tolerance of SysPFMG","HYPRE_SStructSysPFMGSetTol",ex->tol,&ex->tol,NULL);
1974:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(ex->ss_solver,ex->tol));
1975:   PetscOptionsEList("-pc_syspfmg_relax_type","Relax type for the up and down cycles","HYPRE_SStructSysPFMGSetRelaxType",SysPFMGRelaxType,4,SysPFMGRelaxType[ex->relax_type],&ex->relax_type,NULL);
1976:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetRelaxType,(ex->ss_solver, ex->relax_type));
1977:   PetscOptionsTail();
1978:   return(0);
1979: }

1983: PetscErrorCode PCApply_SysPFMG(PC pc,Vec x,Vec y)
1984: {
1985:   PetscErrorCode    ierr;
1986:   PC_SysPFMG        *ex = (PC_SysPFMG*) pc->data;
1987:   PetscScalar       *yy;
1988:   const PetscScalar *xx;
1989:   PetscInt          ilower[3],iupper[3];
1990:   Mat_HYPRESStruct  *mx     = (Mat_HYPRESStruct*)(pc->pmat->data);
1991:   PetscInt          ordering= mx->dofs_order;
1992:   PetscInt          nvars   = mx->nvars;
1993:   PetscInt          part    = 0;
1994:   PetscInt          size;
1995:   PetscInt          i;

1998:   PetscCitationsRegister(hypreCitation,&cite);
1999:   DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);
2000:   iupper[0] += ilower[0] - 1;
2001:   iupper[1] += ilower[1] - 1;
2002:   iupper[2] += ilower[2] - 1;

2004:   size = 1;
2005:   for (i= 0; i< 3; i++) size *= (iupper[i]-ilower[i]+1);

2007:   /* copy x values over to hypre for variable ordering */
2008:   if (ordering) {
2009:     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
2010:     VecGetArrayRead(x,&xx);
2011:     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,(PetscScalar*)xx+(size*i)));
2012:     VecRestoreArrayRead(x,&xx);
2013:     PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
2014:     PetscStackCallStandard(HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
2015:     PetscStackCallStandard(HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));

2017:     /* copy solution values back to PETSc */
2018:     VecGetArray(y,&yy);
2019:     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,yy+(size*i)));
2020:     VecRestoreArray(y,&yy);
2021:   } else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
2022:     PetscScalar *z;
2023:     PetscInt    j, k;

2025:     PetscMalloc1(nvars*size,&z);
2026:     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
2027:     VecGetArrayRead(x,&xx);

2029:     /* transform nodal to hypre's variable ordering for sys_pfmg */
2030:     for (i= 0; i< size; i++) {
2031:       k= i*nvars;
2032:       for (j= 0; j< nvars; j++) z[j*size+i]= xx[k+j];
2033:     }
2034:     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,z+(size*i)));
2035:     VecRestoreArrayRead(x,&xx);
2036:     PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
2037:     PetscStackCallStandard(HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));

2039:     /* copy solution values back to PETSc */
2040:     VecGetArray(y,&yy);
2041:     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,z+(size*i)));
2042:     /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
2043:     for (i= 0; i< size; i++) {
2044:       k= i*nvars;
2045:       for (j= 0; j< nvars; j++) yy[k+j]= z[j*size+i];
2046:     }
2047:     VecRestoreArray(y,&yy);
2048:     PetscFree(z);
2049:   }
2050:   return(0);
2051: }

2055: 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)
2056: {
2057:   PC_SysPFMG     *jac = (PC_SysPFMG*)pc->data;
2059:   PetscInt       oits;

2062:   PetscCitationsRegister(hypreCitation,&cite);
2063:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,its*jac->its));
2064:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,rtol));
2065:   PCApply_SysPFMG(pc,b,y);
2066:   PetscStackCallStandard(HYPRE_SStructSysPFMGGetNumIterations,(jac->ss_solver,(HYPRE_Int *)&oits));
2067:   *outits = oits;
2068:   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
2069:   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
2070:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,jac->tol));
2071:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,jac->its));
2072:   return(0);
2073: }


2078: PetscErrorCode PCSetUp_SysPFMG(PC pc)
2079: {
2080:   PetscErrorCode   ierr;
2081:   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
2082:   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct*)(pc->pmat->data);
2083:   PetscBool        flg;

2086:   PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESSTRUCT,&flg);
2087:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESSTRUCT with this preconditioner");

2089:   /* create the hypre sstruct solver object and set its information */
2090:   if (ex->ss_solver) PetscStackCallStandard(HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));
2091:   PetscStackCallStandard(HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
2092:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetZeroGuess,(ex->ss_solver));
2093:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetup,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
2094:   return(0);
2095: }


2098: /*MC
2099:      PCSysPFMG - the hypre SysPFMG multigrid solver

2101:    Level: advanced

2103:    Options Database:
2104: + -pc_syspfmg_its <its> number of iterations of SysPFMG to use as preconditioner
2105: . -pc_syspfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
2106: . -pc_syspfmg_num_post_relax <steps> number of smoothing steps after coarse grid
2107: . -pc_syspfmg_tol <tol> tolerance of SysPFMG
2108: . -pc_syspfmg_relax_type -relaxation type for the up and down cycles, one of Weighted-Jacobi,Red/Black-Gauss-Seidel

2110:    Notes:  This is for CELL-centered descretizations

2112:            This must be used with the MATHYPRESSTRUCT matrix type.
2113:            This is less general than in hypre, it supports only one part, and one block per process defined by a PETSc DMDA.
2114:            Also, only cell-centered variables.

2116: .seealso:  PCMG, MATHYPRESSTRUCT
2117: M*/

2121: PETSC_EXTERN PetscErrorCode PCCreate_SysPFMG(PC pc)
2122: {
2124:   PC_SysPFMG     *ex;

2127:   PetscNew(&ex); \
2128:   pc->data = ex;

2130:   ex->its            = 1;
2131:   ex->tol            = 1.e-8;
2132:   ex->relax_type     = 1;
2133:   ex->num_pre_relax  = 1;
2134:   ex->num_post_relax = 1;

2136:   pc->ops->setfromoptions  = PCSetFromOptions_SysPFMG;
2137:   pc->ops->view            = PCView_SysPFMG;
2138:   pc->ops->destroy         = PCDestroy_SysPFMG;
2139:   pc->ops->apply           = PCApply_SysPFMG;
2140:   pc->ops->applyrichardson = PCApplyRichardson_SysPFMG;
2141:   pc->ops->setup           = PCSetUp_SysPFMG;

2143:   MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ex->hcomm));
2144:   PetscStackCallStandard(HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
2145:   return(0);
2146: }