Actual source code: hypre.c

petsc-master 2017-02-27
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>
  9: /* this include is needed ONLY to allow access to the private data inside the Mat object specific to hypre */
 10:  #include <petsc/private/matimpl.h>
 11:  #include <../src/vec/vec/impls/hypre/vhyp.h>
 12:  #include <../src/mat/impls/hypre/mhypre.h>
 13:  #include <../src/dm/impls/da/hypre/mhyp.h>
 14: #include <_hypre_parcsr_ls.h>

 16: static PetscBool cite = PETSC_FALSE;
 17: 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";

 19: /*
 20:    Private context (data structure) for the  preconditioner.
 21: */
 22: typedef struct {
 23:   HYPRE_Solver   hsolver;
 24:   Mat            hpmat; /* MatHYPRE */

 26:   HYPRE_Int (*destroy)(HYPRE_Solver);
 27:   HYPRE_Int (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
 28:   HYPRE_Int (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);

 30:   MPI_Comm comm_hypre;
 31:   char     *hypre_type;

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

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

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

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

 53:   /* options for BoomerAMG */
 54:   PetscInt  cycletype;
 55:   PetscInt  maxlevels;
 56:   double    strongthreshold;
 57:   double    maxrowsum;
 58:   PetscInt  gridsweeps[3];
 59:   PetscInt  coarsentype;
 60:   PetscInt  measuretype;
 61:   PetscInt  smoothtype;
 62:   PetscInt  smoothnumlevels;
 63:   PetscInt  eu_level;   /* Number of levels for ILU(k) in Euclid */
 64:   double    eu_droptolerance; /* Drop tolerance for ILU(k) in Euclid */
 65:   PetscInt  eu_bj;      /* Defines use of Block Jacobi ILU in Euclid */
 66:   PetscInt  relaxtype[3];
 67:   double    relaxweight;
 68:   double    outerrelaxweight;
 69:   PetscInt  relaxorder;
 70:   double    truncfactor;
 71:   PetscBool applyrichardson;
 72:   PetscInt  pmax;
 73:   PetscInt  interptype;
 74:   PetscInt  agg_nl;
 75:   PetscInt  agg_num_paths;
 76:   PetscInt  nodal_coarsen;
 77:   PetscBool nodal_relax;
 78:   PetscInt  nodal_relax_levels;

 80:   PetscInt  nodal_coarsening;
 81:   PetscInt  vec_interp_variant;
 82:   HYPRE_IJVector  *hmnull;
 83:   HYPRE_ParVector *phmnull;  /* near null space passed to hypre */
 84:   PetscInt        n_hmnull;
 85:   Vec             hmnull_constant;
 86:   PetscScalar     **hmnull_hypre_data_array;   /* this is the space in hmnull that was allocated by hypre, it is restored to hypre just before freeing the phmnull vectors */

 88:   /* options for AS (Auxiliary Space preconditioners) */
 89:   PetscInt  as_print;
 90:   PetscInt  as_max_iter;
 91:   PetscReal as_tol;
 92:   PetscInt  as_relax_type;
 93:   PetscInt  as_relax_times;
 94:   PetscReal as_relax_weight;
 95:   PetscReal as_omega;
 96:   PetscInt  as_amg_alpha_opts[5]; /* AMG coarsen type, agg_levels, relax_type, interp_type, Pmax for vector Poisson (AMS) or Curl problem (ADS) */
 97:   PetscReal as_amg_alpha_theta;   /* AMG strength for vector Poisson (AMS) or Curl problem (ADS) */
 98:   PetscInt  as_amg_beta_opts[5];  /* AMG coarsen type, agg_levels, relax_type, interp_type, Pmax for scalar Poisson (AMS) or vector Poisson (ADS) */
 99:   PetscReal as_amg_beta_theta;    /* AMG strength for scalar Poisson (AMS) or vector Poisson (ADS)  */
100:   PetscInt  ams_cycle_type;
101:   PetscInt  ads_cycle_type;

103:   /* additional data */
104:   Mat G;             /* MatHYPRE */
105:   Mat C;             /* MatHYPRE */
106:   Mat alpha_Poisson; /* MatHYPRE */
107:   Mat beta_Poisson;  /* MatHYPRE */

109:   /* extra information for AMS */
110:   PetscInt       dim; /* geometrical dimension */
111:   HYPRE_IJVector coords[3];
112:   HYPRE_IJVector constants[3];
113:   PetscBool      ams_beta_is_zero;
114:   PetscBool      ams_beta_is_zero_part;
115:   PetscInt       ams_proj_freq;
116: } PC_HYPRE;

118: PetscErrorCode PCHYPREGetSolver(PC pc,HYPRE_Solver *hsolver)
119: {
120:   PC_HYPRE *jac = (PC_HYPRE*)pc->data;

123:   *hsolver = jac->hsolver;
124:   return(0);
125: }

127: static PetscErrorCode PCSetUp_HYPRE(PC pc)
128: {
129:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
130:   Mat_HYPRE          *hjac;
131:   HYPRE_ParCSRMatrix hmat;
132:   HYPRE_ParVector    bv,xv;
133:   PetscBool          ishypre;
134:   PetscErrorCode     ierr;

137:   if (!jac->hypre_type) {
138:     PCHYPRESetType(pc,"boomeramg");
139:   }

141:   PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRE,&ishypre);
142:   if (!ishypre) {
143:     MatReuse reuse;
144:     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
145:     else reuse = MAT_INITIAL_MATRIX;
146:     MatConvert(pc->pmat,MATHYPRE,reuse,&jac->hpmat);
147:   } else {
148:     PetscObjectReference((PetscObject)pc->pmat);
149:     MatDestroy(&jac->hpmat);
150:     jac->hpmat = pc->pmat;
151:   }
152:   hjac = (Mat_HYPRE*)(jac->hpmat->data);

154:   /* special case for BoomerAMG */
155:   if (jac->setup == HYPRE_BoomerAMGSetup) {
156:     MatNullSpace    mnull;
157:     PetscBool       has_const;
158:     PetscInt        bs,nvec,i;
159:     const Vec       *vecs;
160:     PetscScalar     *petscvecarray;
161: 
162:     MatGetBlockSize(pc->pmat,&bs);
163:     if (bs > 1) PetscStackCallStandard(HYPRE_BoomerAMGSetNumFunctions,(jac->hsolver,bs));
164:     MatGetNearNullSpace(pc->mat, &mnull);
165:     if (mnull) {
166:       MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);
167:       PetscMalloc1(nvec+1,&jac->hmnull);
168:       PetscMalloc1(nvec+1,&jac->hmnull_hypre_data_array);
169:       PetscMalloc1(nvec+1,&jac->phmnull);
170:       for (i=0; i<nvec; i++) {
171:         VecHYPRE_IJVectorCreate(vecs[i],&jac->hmnull[i]);
172:         VecGetArrayRead(vecs[i],(const PetscScalar **)&petscvecarray);
173:         VecHYPRE_ParVectorReplacePointer(jac->hmnull[i],petscvecarray,jac->hmnull_hypre_data_array[i]);
174:         VecRestoreArrayRead(vecs[i],(const PetscScalar **)&petscvecarray);
175:         PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->hmnull[i],(void**)&jac->phmnull[i]));
176:       }
177:       if (has_const) {
178:         MatCreateVecs(pc->pmat,&jac->hmnull_constant,NULL);
179:         VecSet(jac->hmnull_constant,1);
180:         VecNormalize(jac->hmnull_constant,NULL);
181:         VecHYPRE_IJVectorCreate(jac->hmnull_constant,&jac->hmnull[nvec]);
182:         VecGetArrayRead(jac->hmnull_constant,(const PetscScalar **)&petscvecarray);
183:         VecHYPRE_ParVectorReplacePointer(jac->hmnull[nvec],petscvecarray,jac->hmnull_hypre_data_array[nvec]);
184:         VecRestoreArrayRead(jac->hmnull_constant,(const PetscScalar **)&petscvecarray);
185:         PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->hmnull[nvec],(void**)&jac->phmnull[nvec]));
186:         nvec++;
187:       }
188:       PetscStackCallStandard(HYPRE_BoomerAMGSetInterpVectors,(jac->hsolver,nvec,jac->phmnull));
189:       jac->n_hmnull = nvec;
190:     }
191:   }

193:   /* special case for AMS */
194:   if (jac->setup == HYPRE_AMSSetup) {
195:     Mat_HYPRE          *hm;
196:     HYPRE_ParCSRMatrix parcsr;
197:     if (!jac->coords[0] && !jac->constants[0]) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE AMS preconditioner needs either the coordinate vectors via PCSetCoordinates() or the edge constant vectors via PCHYPRESetEdgeConstantVectors()");
198:     if (jac->dim) {
199:       PetscStackCallStandard(HYPRE_AMSSetDimension,(jac->hsolver,jac->dim));
200:     }
201:     if (jac->constants[0]) {
202:       HYPRE_ParVector ozz,zoz,zzo = NULL;
203:       PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->constants[0],(void**)(&ozz)));
204:       PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->constants[1],(void**)(&zoz)));
205:       if (jac->constants[2]) {
206:         PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->constants[2],(void**)(&zzo)));
207:       }
208:       PetscStackCallStandard(HYPRE_AMSSetEdgeConstantVectors,(jac->hsolver,ozz,zoz,zzo));
209:     }
210:     if (jac->coords[0]) {
211:       HYPRE_ParVector coords[3];
212:       coords[0] = NULL;
213:       coords[1] = NULL;
214:       coords[2] = NULL;
215:       if (jac->coords[0]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[0],(void**)(&coords[0])));
216:       if (jac->coords[1]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[1],(void**)(&coords[1])));
217:       if (jac->coords[2]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[2],(void**)(&coords[2])));
218:       PetscStackCallStandard(HYPRE_AMSSetCoordinateVectors,(jac->hsolver,coords[0],coords[1],coords[2]));
219:     }
220:     if (!jac->G) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE AMS preconditioner needs the discrete gradient operator via PCHYPRESetDiscreteGradient");
221:     hm = (Mat_HYPRE*)(jac->G->data);
222:     PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&parcsr)));
223:     PetscStackCallStandard(HYPRE_AMSSetDiscreteGradient,(jac->hsolver,parcsr));
224:     if (jac->alpha_Poisson) {
225:       hm = (Mat_HYPRE*)(jac->alpha_Poisson->data);
226:       PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&parcsr)));
227:       PetscStackCallStandard(HYPRE_AMSSetAlphaPoissonMatrix,(jac->hsolver,parcsr));
228:     }
229:     if (jac->ams_beta_is_zero) {
230:       PetscStackCallStandard(HYPRE_AMSSetBetaPoissonMatrix,(jac->hsolver,NULL));
231:     } else if (jac->beta_Poisson) {
232:       hm = (Mat_HYPRE*)(jac->beta_Poisson->data);
233:       PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&parcsr)));
234:       PetscStackCallStandard(HYPRE_AMSSetBetaPoissonMatrix,(jac->hsolver,parcsr));
235:     }
236:   }
237:   /* special case for ADS */
238:   if (jac->setup == HYPRE_ADSSetup) {
239:     Mat_HYPRE          *hm;
240:     HYPRE_ParCSRMatrix parcsr;
241:     if (!jac->coords[0]) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE ADS preconditioner needs the coordinate vectors via PCSetCoordinates()");
242:     else if (!jac->coords[1] || !jac->coords[2]) 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");
243:     if (!jac->G) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE ADS preconditioner needs the discrete gradient operator via PCHYPRESetDiscreteGradient");
244:     if (!jac->C) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE ADS preconditioner needs the discrete curl operator via PCHYPRESetDiscreteGradient");
245:     if (jac->coords[0]) {
246:       HYPRE_ParVector coords[3];
247:       coords[0] = NULL;
248:       coords[1] = NULL;
249:       coords[2] = NULL;
250:       if (jac->coords[0]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[0],(void**)(&coords[0])));
251:       if (jac->coords[1]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[1],(void**)(&coords[1])));
252:       if (jac->coords[2]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[2],(void**)(&coords[2])));
253:       PetscStackCallStandard(HYPRE_ADSSetCoordinateVectors,(jac->hsolver,coords[0],coords[1],coords[2]));
254:     }
255:     hm = (Mat_HYPRE*)(jac->G->data);
256:     PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&parcsr)));
257:     PetscStackCallStandard(HYPRE_ADSSetDiscreteGradient,(jac->hsolver,parcsr));
258:     hm = (Mat_HYPRE*)(jac->C->data);
259:     PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&parcsr)));
260:     PetscStackCallStandard(HYPRE_ADSSetDiscreteCurl,(jac->hsolver,parcsr));
261:   }
262:   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hjac->ij,(void**)&hmat));
263:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hjac->b,(void**)&bv));
264:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hjac->x,(void**)&xv));
265:   PetscStackCall("HYPRE_SetupXXX",(*jac->setup)(jac->hsolver,hmat,bv,xv););
266:   return(0);
267: }

269: static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
270: {
271:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
272:   Mat_HYPRE          *hjac = (Mat_HYPRE*)(jac->hpmat->data);
273:   PetscErrorCode     ierr;
274:   HYPRE_ParCSRMatrix hmat;
275:   PetscScalar        *xv;
276:   const PetscScalar  *bv,*sbv;
277:   HYPRE_ParVector    jbv,jxv;
278:   PetscScalar        *sxv;
279:   PetscInt           hierr;

282:   PetscCitationsRegister(hypreCitation,&cite);
283:   if (!jac->applyrichardson) {VecSet(x,0.0);}
284:   VecGetArrayRead(b,&bv);
285:   VecGetArray(x,&xv);
286:   VecHYPRE_ParVectorReplacePointer(hjac->b,(PetscScalar*)bv,sbv);
287:   VecHYPRE_ParVectorReplacePointer(hjac->x,xv,sxv);

289:   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hjac->ij,(void**)&hmat));
290:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hjac->b,(void**)&jbv));
291:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hjac->x,(void**)&jxv));
292:   PetscStackCall("Hypre solve",h(*jac->solve)(jac->hsolver,hmat,jbv,jxv);
293:                                if (hierr && hierr != HYPRE_ERROR_CONV) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
294:                                if (hierr) hypre__global_error = 0;);

296:   if (jac->setup == HYPRE_AMSSetup && jac->ams_beta_is_zero_part) {
297:     PetscStackCallStandard(HYPRE_AMSProjectOutGradients,(jac->hsolver,jxv));
298:   }
299:   VecHYPRE_ParVectorReplacePointer(hjac->b,(PetscScalar*)sbv,bv);
300:   VecHYPRE_ParVectorReplacePointer(hjac->x,sxv,xv);
301:   VecRestoreArray(x,&xv);
302:   VecRestoreArrayRead(b,&bv);
303:   return(0);
304: }

306: static PetscErrorCode PCReset_HYPRE(PC pc)
307: {
308:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;

312:   MatDestroy(&jac->hpmat);
313:   MatDestroy(&jac->G);
314:   MatDestroy(&jac->C);
315:   MatDestroy(&jac->alpha_Poisson);
316:   MatDestroy(&jac->beta_Poisson);
317:   if (jac->coords[0]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[0])); jac->coords[0] = NULL;
318:   if (jac->coords[1]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[1])); jac->coords[1] = NULL;
319:   if (jac->coords[2]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[2])); jac->coords[2] = NULL;
320:   if (jac->constants[0]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[0])); jac->constants[0] = NULL;
321:   if (jac->constants[1]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[1])); jac->constants[1] = NULL;
322:   if (jac->constants[2]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[2])); jac->constants[2] = NULL;
323:   if (jac->n_hmnull && jac->hmnull) {
324:     PetscInt                 i;
325:     PETSC_UNUSED PetscScalar *petscvecarray;

327:     for (i=0; i<jac->n_hmnull; i++) {
328:       VecHYPRE_ParVectorReplacePointer(jac->hmnull[i],jac->hmnull_hypre_data_array[i],petscvecarray);
329:       PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->hmnull[i]));
330:     }
331:     PetscFree(jac->hmnull);
332:     PetscFree(jac->hmnull_hypre_data_array);
333:     PetscFree(jac->phmnull);
334:     VecDestroy(&jac->hmnull_constant);
335:   }
336:   jac->ams_beta_is_zero = PETSC_FALSE;
337:   jac->dim = 0;
338:   return(0);
339: }

341: static PetscErrorCode PCDestroy_HYPRE(PC pc)
342: {
343:   PC_HYPRE                 *jac = (PC_HYPRE*)pc->data;
344:   PetscErrorCode           ierr;

347:   PCReset_HYPRE(pc);
348:   if (jac->destroy) PetscStackCall("HYPRE_DestroyXXX",(*jac->destroy)(jac->hsolver););
349:   PetscFree(jac->hypre_type);
350:   if (jac->comm_hypre != MPI_COMM_NULL) { MPI_Comm_free(&(jac->comm_hypre));}
351:   PetscFree(pc->data);

353:   PetscObjectChangeTypeName((PetscObject)pc,0);
354:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetType_C",NULL);
355:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPREGetType_C",NULL);
356:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetCoordinates_C",NULL);
357:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteGradient_C",NULL);
358:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteCurl_C",NULL);
359:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetConstantEdgeVectors_C",NULL);
360:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetPoissonMatrix_C",NULL);
361:   return(0);
362: }

364: /* --------------------------------------------------------------------------------------------*/
365: static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PetscOptionItems *PetscOptionsObject,PC pc)
366: {
367:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
369:   PetscBool      flag;

372:   PetscOptionsHead(PetscOptionsObject,"HYPRE Pilut Options");
373:   PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);
374:   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetMaxIter,(jac->hsolver,jac->maxiter));
375:   PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);
376:   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetDropTolerance,(jac->hsolver,jac->tol));
377:   PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);
378:   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetFactorRowSize,(jac->hsolver,jac->factorrowsize));
379:   PetscOptionsTail();
380:   return(0);
381: }

383: static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
384: {
385:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
387:   PetscBool      iascii;

390:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
391:   if (iascii) {
392:     PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");
393:     if (jac->maxiter != PETSC_DEFAULT) {
394:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);
395:     } else {
396:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");
397:     }
398:     if (jac->tol != PETSC_DEFAULT) {
399:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %g\n",(double)jac->tol);
400:     } else {
401:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");
402:     }
403:     if (jac->factorrowsize != PETSC_DEFAULT) {
404:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);
405:     } else {
406:       PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");
407:     }
408:   }
409:   return(0);
410: }

412: /* --------------------------------------------------------------------------------------------*/

414: static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
415: {
416:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
417:   Mat_HYPRE          *hjac = (Mat_HYPRE*)(jac->hpmat->data);
418:   PetscErrorCode     ierr;
419:   HYPRE_ParCSRMatrix hmat;
420:   PetscScalar        *xv;
421:   const PetscScalar  *bv;
422:   HYPRE_ParVector    jbv,jxv;
423:   PetscScalar        *sbv,*sxv;
424:   PetscInt           hierr;

427:   PetscCitationsRegister(hypreCitation,&cite);
428:   VecSet(x,0.0);
429:   VecGetArrayRead(b,&bv);
430:   VecGetArray(x,&xv);
431:   VecHYPRE_ParVectorReplacePointer(hjac->b,(PetscScalar*)bv,sbv);
432:   VecHYPRE_ParVectorReplacePointer(hjac->x,xv,sxv);

434:   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hjac->ij,(void**)&hmat));
435:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hjac->b,(void**)&jbv));
436:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hjac->x,(void**)&jxv));

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

443:   VecHYPRE_ParVectorReplacePointer(hjac->b,sbv,bv);
444:   VecHYPRE_ParVectorReplacePointer(hjac->x,sxv,xv);
445:   VecRestoreArray(x,&xv);
446:   VecRestoreArrayRead(b,&bv);
447:   return(0);
448: }

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

453: static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
454: static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"};
455: static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
456: /* The following corresponds to HYPRE_BoomerAMGSetRelaxType which has many missing numbers in the enum */
457: static const char *HYPREBoomerAMGSmoothType[]   = {"Schwarz-smoothers","Pilut","ParaSails","Euclid"};
458: static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","seqboundary-Gauss-Seidel","SOR/Jacobi","backward-SOR/Jacobi",
459:                                                   "" /* [5] hybrid chaotic Gauss-Seidel (works only with OpenMP) */,"symmetric-SOR/Jacobi",
460:                                                   "" /* 7 */,"l1scaled-SOR/Jacobi","Gaussian-elimination",
461:                                                   "" /* 10 */, "" /* 11 */, "" /* 12 */, "" /* 13 */, "" /* 14 */,
462:                                                   "CG" /* non-stationary */,"Chebyshev","FCF-Jacobi","l1scaled-Jacobi"};
463: static const char *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i",
464:                                                   "ext+i-cc", "standard", "standard-wts", "block", "block-wtd", "FF", "FF1"};
465: static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PetscOptionItems *PetscOptionsObject,PC pc)
466: {
467:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
469:   PetscInt       n,indx,level;
470:   PetscBool      flg, tmp_truth;
471:   double         tmpdbl, twodbl[2];

474:   PetscOptionsHead(PetscOptionsObject,"HYPRE BoomerAMG Options");
475:   PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType+1,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);
476:   if (flg) {
477:     jac->cycletype = indx+1;
478:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
479:   }
480:   PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);
481:   if (flg) {
482:     if (jac->maxlevels < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
483:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
484:   }
485:   PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);
486:   if (flg) {
487:     if (jac->maxiter < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
488:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
489:   }
490:   PetscOptionsScalar("-pc_hypre_boomeramg_tol","Convergence tolerance PER hypre call (0.0 = use a fixed number of iterations)","None",jac->tol,&jac->tol,&flg);
491:   if (flg) {
492:     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);
493:     PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
494:   }

496:   PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);
497:   if (flg) {
498:     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);
499:     PetscStackCallStandard(HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
500:   }

502:   PetscOptionsInt("-pc_hypre_boomeramg_P_max","Max elements per row for interpolation operator (0=unlimited)","None",jac->pmax,&jac->pmax,&flg);
503:   if (flg) {
504:     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);
505:     PetscStackCallStandard(HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
506:   }

508:   PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);
509:   if (flg) {
510:     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);

512:     PetscStackCallStandard(HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
513:   }


516:   PetscOptionsInt("-pc_hypre_boomeramg_agg_num_paths","Number of paths for aggressive coarsening","None",jac->agg_num_paths,&jac->agg_num_paths,&flg);
517:   if (flg) {
518:     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);

520:     PetscStackCallStandard(HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
521:   }


524:   PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);
525:   if (flg) {
526:     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);
527:     PetscStackCallStandard(HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
528:   }
529:   PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);
530:   if (flg) {
531:     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);
532:     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);
533:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
534:   }

536:   /* Grid sweeps */
537:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_all","Number of sweeps for the up and down grid levels","None",jac->gridsweeps[0],&indx,&flg);
538:   if (flg) {
539:     PetscStackCallStandard(HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver,indx));
540:     /* modify the jac structure so we can view the updated options with PC_View */
541:     jac->gridsweeps[0] = indx;
542:     jac->gridsweeps[1] = indx;
543:     /*defaults coarse to 1 */
544:     jac->gridsweeps[2] = 1;
545:   }

547:   PetscOptionsInt("-pc_hypre_boomeramg_nodal_coarsen","Use a nodal based coarsening 1-6","HYPRE_BoomerAMGSetNodal",jac->nodal_coarsening, &jac->nodal_coarsening,&flg);
548:   if (flg) {
549:     PetscStackCallStandard(HYPRE_BoomerAMGSetNodal,(jac->hsolver,jac->nodal_coarsening));
550:   }

552:   PetscOptionsInt("-pc_hypre_boomeramg_vec_interp_variant","Variant of algorithm 1-3","HYPRE_BoomerAMGSetInterpVecVariant",jac->vec_interp_variant, &jac->vec_interp_variant,&flg);
553:   if (flg) {
554:     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpVecVariant,(jac->hsolver,jac->vec_interp_variant));
555:   }

557:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx,&flg);
558:   if (flg) {
559:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 1));
560:     jac->gridsweeps[0] = indx;
561:   }
562:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);
563:   if (flg) {
564:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 2));
565:     jac->gridsweeps[1] = indx;
566:   }
567:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);
568:   if (flg) {
569:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 3));
570:     jac->gridsweeps[2] = indx;
571:   }

573:   /* Smooth type */
574:   PetscOptionsEList("-pc_hypre_boomeramg_smooth_type","Enable more complex smoothers","None",HYPREBoomerAMGSmoothType,ALEN(HYPREBoomerAMGSmoothType),HYPREBoomerAMGSmoothType[0],&indx,&flg);
575:   if (flg) {
576:     jac->smoothtype = indx;
577:     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothType,(jac->hsolver,indx+6));
578:     jac->smoothnumlevels = 25;
579:     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothNumLevels,(jac->hsolver,25));
580:   }

582:   /* Number of smoothing levels */
583:   PetscOptionsInt("-pc_hypre_boomeramg_smooth_num_levels","Number of levels on which more complex smoothers are used","None",25,&indx,&flg);
584:   if (flg && (jac->smoothtype != -1)) {
585:     jac->smoothnumlevels = indx;
586:     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothNumLevels,(jac->hsolver,indx));
587:   }

589:   /* Number of levels for ILU(k) for Euclid */
590:   PetscOptionsInt("-pc_hypre_boomeramg_eu_level","Number of levels for ILU(k) in Euclid smoother","None",0,&indx,&flg);
591:   if (flg && (jac->smoothtype == 3)) {
592:     jac->eu_level = indx;
593:     PetscStackCallStandard(HYPRE_BoomerAMGSetEuLevel,(jac->hsolver,indx));
594:   }

596:   /* Filter for ILU(k) for Euclid */
597:   double droptolerance;
598:   PetscOptionsScalar("-pc_hypre_boomeramg_eu_droptolerance","Drop tolerance for ILU(k) in Euclid smoother","None",0,&droptolerance,&flg);
599:   if (flg && (jac->smoothtype == 3)) {
600:     jac->eu_droptolerance = droptolerance;
601:     PetscStackCallStandard(HYPRE_BoomerAMGSetEuLevel,(jac->hsolver,droptolerance));
602:   }

604:   /* Use Block Jacobi ILUT for Euclid */
605:   PetscOptionsBool("-pc_hypre_boomeramg_eu_bj", "Use Block Jacobi for ILU in Euclid smoother?", "None", PETSC_FALSE, &tmp_truth, &flg);
606:   if (flg && (jac->smoothtype == 3)) {
607:     jac->eu_bj = tmp_truth;
608:     PetscStackCallStandard(HYPRE_BoomerAMGSetEuBJ,(jac->hsolver,jac->eu_bj));
609:   }

611:   /* Relax type */
612:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_all","Relax type for the up and down cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);
613:   if (flg) {
614:     jac->relaxtype[0] = jac->relaxtype[1]  = indx;
615:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, indx));
616:     /* by default, coarse type set to 9 */
617:     jac->relaxtype[2] = 9;

619:   }
620:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);
621:   if (flg) {
622:     jac->relaxtype[0] = indx;
623:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 1));
624:   }
625:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);
626:   if (flg) {
627:     jac->relaxtype[1] = indx;
628:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 2));
629:   }
630:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[9],&indx,&flg);
631:   if (flg) {
632:     jac->relaxtype[2] = indx;
633:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 3));
634:   }

636:   /* Relaxation Weight */
637:   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);
638:   if (flg) {
639:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxWt,(jac->hsolver,tmpdbl));
640:     jac->relaxweight = tmpdbl;
641:   }

643:   n         = 2;
644:   twodbl[0] = twodbl[1] = 1.0;
645:   PetscOptionsRealArray("-pc_hypre_boomeramg_relax_weight_level","Set the relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);
646:   if (flg) {
647:     if (n == 2) {
648:       indx =  (int)PetscAbsReal(twodbl[1]);
649:       PetscStackCallStandard(HYPRE_BoomerAMGSetLevelRelaxWt,(jac->hsolver,twodbl[0],indx));
650:     } 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);
651:   }

653:   /* Outer relaxation Weight */
654:   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);
655:   if (flg) {
656:     PetscStackCallStandard(HYPRE_BoomerAMGSetOuterWt,(jac->hsolver, tmpdbl));
657:     jac->outerrelaxweight = tmpdbl;
658:   }

660:   n         = 2;
661:   twodbl[0] = twodbl[1] = 1.0;
662:   PetscOptionsRealArray("-pc_hypre_boomeramg_outer_relax_weight_level","Set the outer relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);
663:   if (flg) {
664:     if (n == 2) {
665:       indx =  (int)PetscAbsReal(twodbl[1]);
666:       PetscStackCallStandard(HYPRE_BoomerAMGSetLevelOuterWt,(jac->hsolver, twodbl[0], indx));
667:     } 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);
668:   }

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

673:   if (flg && tmp_truth) {
674:     jac->relaxorder = 0;
675:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
676:   }
677:   PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,ALEN(HYPREBoomerAMGMeasureType),HYPREBoomerAMGMeasureType[0],&indx,&flg);
678:   if (flg) {
679:     jac->measuretype = indx;
680:     PetscStackCallStandard(HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
681:   }
682:   /* update list length 3/07 */
683:   PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,ALEN(HYPREBoomerAMGCoarsenType),HYPREBoomerAMGCoarsenType[6],&indx,&flg);
684:   if (flg) {
685:     jac->coarsentype = indx;
686:     PetscStackCallStandard(HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
687:   }

689:   /* new 3/07 */
690:   PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,ALEN(HYPREBoomerAMGInterpType),HYPREBoomerAMGInterpType[0],&indx,&flg);
691:   if (flg) {
692:     jac->interptype = indx;
693:     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
694:   }

696:   PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);
697:   if (flg) {
698:     level = 3;
699:     PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,NULL);

701:     jac->printstatistics = PETSC_TRUE;
702:     PetscStackCallStandard(HYPRE_BoomerAMGSetPrintLevel,(jac->hsolver,level));
703:   }

705:   PetscOptionsName("-pc_hypre_boomeramg_print_debug","Print debug information","None",&flg);
706:   if (flg) {
707:     level = 3;
708:     PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,NULL);

710:     jac->printstatistics = PETSC_TRUE;
711:     PetscStackCallStandard(HYPRE_BoomerAMGSetDebugFlag,(jac->hsolver,level));
712:   }

714:   PetscOptionsBool("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);
715:   if (flg && tmp_truth) {
716:     PetscInt tmp_int;
717:     PetscOptionsInt("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);
718:     if (flg) jac->nodal_relax_levels = tmp_int;
719:     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothType,(jac->hsolver,6));
720:     PetscStackCallStandard(HYPRE_BoomerAMGSetDomainType,(jac->hsolver,1));
721:     PetscStackCallStandard(HYPRE_BoomerAMGSetOverlap,(jac->hsolver,0));
722:     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothNumLevels,(jac->hsolver,jac->nodal_relax_levels));
723:   }

725:   PetscOptionsTail();
726:   return(0);
727: }

729: 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)
730: {
731:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
733:   PetscInt       oits;

736:   PetscCitationsRegister(hypreCitation,&cite);
737:   PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,its*jac->maxiter));
738:   PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,rtol));
739:   jac->applyrichardson = PETSC_TRUE;
740:   PCApply_HYPRE(pc,b,y);
741:   jac->applyrichardson = PETSC_FALSE;
742:   PetscStackCallStandard(HYPRE_BoomerAMGGetNumIterations,(jac->hsolver,(HYPRE_Int *)&oits));
743:   *outits = oits;
744:   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
745:   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
746:   PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
747:   PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
748:   return(0);
749: }


752: static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
753: {
754:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
756:   PetscBool      iascii;

759:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
760:   if (iascii) {
761:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");
762:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);
763:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);
764:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);
765:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance PER hypre call %g\n",(double)jac->tol);
766:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %g\n",(double)jac->strongthreshold);
767:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation truncation factor %g\n",(double)jac->truncfactor);
768:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);
769:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);
770:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);

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

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

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

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

785:     if (jac->relaxorder) {
786:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using CF-relaxation\n");
787:     } else {
788:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using CF-relaxation\n");
789:     }
790:     if (jac->smoothtype!=-1) {
791:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Smooth type          %s\n",HYPREBoomerAMGSmoothType[jac->smoothtype]);
792:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Smooth num levels    %d\n",jac->smoothnumlevels);
793:     } else {
794:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using more complex smoothers.\n");
795:     }
796:     if (jac->smoothtype==3) {
797:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Euclid ILU(k) levels %d\n",jac->eu_level);
798:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Euclid ILU(k) drop tolerance %g\n",jac->eu_droptolerance);
799:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Euclid ILU use Block-Jacobi? %d\n",jac->eu_bj);
800:     }
801:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);
802:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);
803:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation type  %s\n",HYPREBoomerAMGInterpType[jac->interptype]);
804:     if (jac->nodal_coarsening) {
805:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal() %D\n",jac->nodal_coarsening);
806:     }
807:     if (jac->vec_interp_variant) {
808:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: HYPRE_BoomerAMGSetInterpVecVariant() %D\n",jac->vec_interp_variant);
809:     }
810:     if (jac->nodal_relax) {
811:       PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);
812:     }
813:   }
814:   return(0);
815: }

817: /* --------------------------------------------------------------------------------------------*/
818: static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PetscOptionItems *PetscOptionsObject,PC pc)
819: {
820:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
822:   PetscInt       indx;
823:   PetscBool      flag;
824:   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};

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

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

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

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

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

844:   PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,ALEN(symtlist),symtlist[0],&indx,&flag);
845:   if (flag) {
846:     jac->symt = indx;
847:     PetscStackCallStandard(HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
848:   }

850:   PetscOptionsTail();
851:   return(0);
852: }

854: static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
855: {
856:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
858:   PetscBool      iascii;
859:   const char     *symt = 0;;

862:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
863:   if (iascii) {
864:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");
865:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);
866:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %g\n",(double)jac->threshhold);
867:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %g\n",(double)jac->filter);
868:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %g\n",(double)jac->loadbal);
869:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscBools[jac->ruse]);
870:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscBools[jac->logging]);
871:     if (!jac->symt) symt = "nonsymmetric matrix and preconditioner";
872:     else if (jac->symt == 1) symt = "SPD matrix and preconditioner";
873:     else if (jac->symt == 2) symt = "nonsymmetric matrix but SPD preconditioner";
874:     else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
875:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);
876:   }
877:   return(0);
878: }
879: /* --------------------------------------------------------------------------------------------*/
880: static PetscErrorCode PCSetFromOptions_HYPRE_AMS(PetscOptionItems *PetscOptionsObject,PC pc)
881: {
882:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
884:   PetscInt       n;
885:   PetscBool      flag,flag2,flag3,flag4;

888:   PetscOptionsHead(PetscOptionsObject,"HYPRE AMS Options");
889:   PetscOptionsInt("-pc_hypre_ams_print_level","Debugging output level for AMS","None",jac->as_print,&jac->as_print,&flag);
890:   if (flag) PetscStackCallStandard(HYPRE_AMSSetPrintLevel,(jac->hsolver,jac->as_print));
891:   PetscOptionsInt("-pc_hypre_ams_max_iter","Maximum number of AMS multigrid iterations within PCApply","None",jac->as_max_iter,&jac->as_max_iter,&flag);
892:   if (flag) PetscStackCallStandard(HYPRE_AMSSetMaxIter,(jac->hsolver,jac->as_max_iter));
893:   PetscOptionsInt("-pc_hypre_ams_cycle_type","Cycle type for AMS multigrid","None",jac->ams_cycle_type,&jac->ams_cycle_type,&flag);
894:   if (flag) PetscStackCallStandard(HYPRE_AMSSetCycleType,(jac->hsolver,jac->ams_cycle_type));
895:   PetscOptionsReal("-pc_hypre_ams_tol","Error tolerance for AMS multigrid","None",jac->as_tol,&jac->as_tol,&flag);
896:   if (flag) PetscStackCallStandard(HYPRE_AMSSetTol,(jac->hsolver,jac->as_tol));
897:   PetscOptionsInt("-pc_hypre_ams_relax_type","Relaxation type for AMS smoother","None",jac->as_relax_type,&jac->as_relax_type,&flag);
898:   PetscOptionsInt("-pc_hypre_ams_relax_times","Number of relaxation steps for AMS smoother","None",jac->as_relax_times,&jac->as_relax_times,&flag2);
899:   PetscOptionsReal("-pc_hypre_ams_relax_weight","Relaxation weight for AMS smoother","None",jac->as_relax_weight,&jac->as_relax_weight,&flag3);
900:   PetscOptionsReal("-pc_hypre_ams_omega","SSOR coefficient for AMS smoother","None",jac->as_omega,&jac->as_omega,&flag4);
901:   if (flag || flag2 || flag3 || flag4) {
902:     PetscStackCallStandard(HYPRE_AMSSetSmoothingOptions,(jac->hsolver,jac->as_relax_type,
903:                                                                       jac->as_relax_times,
904:                                                                       jac->as_relax_weight,
905:                                                                       jac->as_omega));
906:   }
907:   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);
908:   n = 5;
909:   PetscOptionsIntArray("-pc_hypre_ams_amg_alpha_options","AMG options for vector Poisson","None",jac->as_amg_alpha_opts,&n,&flag2);
910:   if (flag || flag2) {
911:     PetscStackCallStandard(HYPRE_AMSSetAlphaAMGOptions,(jac->hsolver,jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
912:                                                                      jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
913:                                                                      jac->as_amg_alpha_opts[2],       /* AMG relax_type */
914:                                                                      jac->as_amg_alpha_theta,
915:                                                                      jac->as_amg_alpha_opts[3],       /* AMG interp_type */
916:                                                                      jac->as_amg_alpha_opts[4]));     /* AMG Pmax */
917:   }
918:   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);
919:   n = 5;
920:   PetscOptionsIntArray("-pc_hypre_ams_amg_beta_options","AMG options for scalar Poisson solver","None",jac->as_amg_beta_opts,&n,&flag2);
921:   if (flag || flag2) {
922:     PetscStackCallStandard(HYPRE_AMSSetBetaAMGOptions,(jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
923:                                                                     jac->as_amg_beta_opts[1],       /* AMG agg_levels */
924:                                                                     jac->as_amg_beta_opts[2],       /* AMG relax_type */
925:                                                                     jac->as_amg_beta_theta,
926:                                                                     jac->as_amg_beta_opts[3],       /* AMG interp_type */
927:                                                                     jac->as_amg_beta_opts[4]));     /* AMG Pmax */
928:   }
929:   PetscOptionsInt("-pc_hypre_ams_projection_frequency","Frequency at which a projection onto the compatible subspace for problems with zero conductivity regions is performed","None",jac->ams_proj_freq,&jac->ams_proj_freq,&flag);
930:   if (flag) { /* override HYPRE's default only if the options is used */
931:     PetscStackCallStandard(HYPRE_AMSSetProjectionFrequency,(jac->hsolver,jac->ams_proj_freq));
932:   }
933:   PetscOptionsTail();
934:   return(0);
935: }

937: static PetscErrorCode PCView_HYPRE_AMS(PC pc,PetscViewer viewer)
938: {
939:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
941:   PetscBool      iascii;

944:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
945:   if (iascii) {
946:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS preconditioning\n");
947:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: subspace iterations per application %d\n",jac->as_max_iter);
948:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: subspace cycle type %d\n",jac->ams_cycle_type);
949:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: subspace iteration tolerance %g\n",jac->as_tol);
950:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: smoother type %d\n",jac->as_relax_type);
951:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: number of smoothing steps %d\n",jac->as_relax_times);
952:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: smoother weight %g\n",jac->as_relax_weight);
953:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: smoother omega %g\n",jac->as_omega);
954:     if (jac->alpha_Poisson) {
955:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: vector Poisson solver (passed in by user)\n");
956:     } else {
957:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: vector Poisson solver (computed) \n");
958:     }
959:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG coarsening type %d\n",jac->as_amg_alpha_opts[0]);
960:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG levels of aggressive coarsening %d\n",jac->as_amg_alpha_opts[1]);
961:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG relaxation type %d\n",jac->as_amg_alpha_opts[2]);
962:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG interpolation type %d\n",jac->as_amg_alpha_opts[3]);
963:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG max nonzero elements in interpolation rows %d\n",jac->as_amg_alpha_opts[4]);
964:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG strength threshold %g\n",jac->as_amg_alpha_theta);
965:     if (!jac->ams_beta_is_zero) {
966:       if (jac->beta_Poisson) {
967:         PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: scalar Poisson solver (passed in by user)\n");
968:       } else {
969:         PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: scalar Poisson solver (computed) \n");
970:       }
971:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG coarsening type %d\n",jac->as_amg_beta_opts[0]);
972:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG levels of aggressive coarsening %d\n",jac->as_amg_beta_opts[1]);
973:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG relaxation type %d\n",jac->as_amg_beta_opts[2]);
974:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG interpolation type %d\n",jac->as_amg_beta_opts[3]);
975:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG max nonzero elements in interpolation rows %d\n",jac->as_amg_beta_opts[4]);
976:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     boomerAMG strength threshold %g\n",jac->as_amg_beta_theta);
977:       if (jac->ams_beta_is_zero_part) {
978:         PetscViewerASCIIPrintf(viewer,"  HYPRE AMS:     compatible subspace projection frequency %d (-1 HYPRE uses default)\n",jac->ams_proj_freq);
979:       }
980:     } else {
981:       PetscViewerASCIIPrintf(viewer,"  HYPRE AMS: scalar Poisson solver not used (zero-conductivity everywhere) \n");
982:     }
983:   }
984:   return(0);
985: }

987: static PetscErrorCode PCSetFromOptions_HYPRE_ADS(PetscOptionItems *PetscOptionsObject,PC pc)
988: {
989:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
991:   PetscInt       n;
992:   PetscBool      flag,flag2,flag3,flag4;

995:   PetscOptionsHead(PetscOptionsObject,"HYPRE ADS Options");
996:   PetscOptionsInt("-pc_hypre_ads_print_level","Debugging output level for ADS","None",jac->as_print,&jac->as_print,&flag);
997:   if (flag) PetscStackCallStandard(HYPRE_ADSSetPrintLevel,(jac->hsolver,jac->as_print));
998:   PetscOptionsInt("-pc_hypre_ads_max_iter","Maximum number of ADS multigrid iterations within PCApply","None",jac->as_max_iter,&jac->as_max_iter,&flag);
999:   if (flag) PetscStackCallStandard(HYPRE_ADSSetMaxIter,(jac->hsolver,jac->as_max_iter));
1000:   PetscOptionsInt("-pc_hypre_ads_cycle_type","Cycle type for ADS multigrid","None",jac->ads_cycle_type,&jac->ads_cycle_type,&flag);
1001:   if (flag) PetscStackCallStandard(HYPRE_ADSSetCycleType,(jac->hsolver,jac->ads_cycle_type));
1002:   PetscOptionsReal("-pc_hypre_ads_tol","Error tolerance for ADS multigrid","None",jac->as_tol,&jac->as_tol,&flag);
1003:   if (flag) PetscStackCallStandard(HYPRE_ADSSetTol,(jac->hsolver,jac->as_tol));
1004:   PetscOptionsInt("-pc_hypre_ads_relax_type","Relaxation type for ADS smoother","None",jac->as_relax_type,&jac->as_relax_type,&flag);
1005:   PetscOptionsInt("-pc_hypre_ads_relax_times","Number of relaxation steps for ADS smoother","None",jac->as_relax_times,&jac->as_relax_times,&flag2);
1006:   PetscOptionsReal("-pc_hypre_ads_relax_weight","Relaxation weight for ADS smoother","None",jac->as_relax_weight,&jac->as_relax_weight,&flag3);
1007:   PetscOptionsReal("-pc_hypre_ads_omega","SSOR coefficient for ADS smoother","None",jac->as_omega,&jac->as_omega,&flag4);
1008:   if (flag || flag2 || flag3 || flag4) {
1009:     PetscStackCallStandard(HYPRE_ADSSetSmoothingOptions,(jac->hsolver,jac->as_relax_type,
1010:                                                                       jac->as_relax_times,
1011:                                                                       jac->as_relax_weight,
1012:                                                                       jac->as_omega));
1013:   }
1014:   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);
1015:   n = 5;
1016:   PetscOptionsIntArray("-pc_hypre_ads_ams_options","AMG options for AMS solver inside ADS","None",jac->as_amg_alpha_opts,&n,&flag2);
1017:   PetscOptionsInt("-pc_hypre_ads_ams_cycle_type","Cycle type for AMS solver inside ADS","None",jac->ams_cycle_type,&jac->ams_cycle_type,&flag3);
1018:   if (flag || flag2 || flag3) {
1019:     PetscStackCallStandard(HYPRE_ADSSetAMSOptions,(jac->hsolver,jac->ams_cycle_type,             /* AMS cycle type */
1020:                                                                 jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
1021:                                                                 jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
1022:                                                                 jac->as_amg_alpha_opts[2],       /* AMG relax_type */
1023:                                                                 jac->as_amg_alpha_theta,
1024:                                                                 jac->as_amg_alpha_opts[3],       /* AMG interp_type */
1025:                                                                 jac->as_amg_alpha_opts[4]));     /* AMG Pmax */
1026:   }
1027:   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);
1028:   n = 5;
1029:   PetscOptionsIntArray("-pc_hypre_ads_amg_options","AMG options for vector AMG solver inside ADS","None",jac->as_amg_beta_opts,&n,&flag2);
1030:   if (flag || flag2) {
1031:     PetscStackCallStandard(HYPRE_ADSSetAMGOptions,(jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
1032:                                                                 jac->as_amg_beta_opts[1],       /* AMG agg_levels */
1033:                                                                 jac->as_amg_beta_opts[2],       /* AMG relax_type */
1034:                                                                 jac->as_amg_beta_theta,
1035:                                                                 jac->as_amg_beta_opts[3],       /* AMG interp_type */
1036:                                                                 jac->as_amg_beta_opts[4]));     /* AMG Pmax */
1037:   }
1038:   PetscOptionsTail();
1039:   return(0);
1040: }

1042: static PetscErrorCode PCView_HYPRE_ADS(PC pc,PetscViewer viewer)
1043: {
1044:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1046:   PetscBool      iascii;

1049:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1050:   if (iascii) {
1051:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS preconditioning\n");
1052:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: subspace iterations per application %d\n",jac->as_max_iter);
1053:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: subspace cycle type %d\n",jac->ads_cycle_type);
1054:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: subspace iteration tolerance %g\n",jac->as_tol);
1055:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: smoother type %d\n",jac->as_relax_type);
1056:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: number of smoothing steps %d\n",jac->as_relax_times);
1057:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: smoother weight %g\n",jac->as_relax_weight);
1058:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: smoother omega %g\n",jac->as_omega);
1059:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: AMS solver\n");
1060:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     subspace cycle type %d\n",jac->ams_cycle_type);
1061:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG coarsening type %d\n",jac->as_amg_alpha_opts[0]);
1062:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG levels of aggressive coarsening %d\n",jac->as_amg_alpha_opts[1]);
1063:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG relaxation type %d\n",jac->as_amg_alpha_opts[2]);
1064:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG interpolation type %d\n",jac->as_amg_alpha_opts[3]);
1065:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG max nonzero elements in interpolation rows %d\n",jac->as_amg_alpha_opts[4]);
1066:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG strength threshold %g\n",jac->as_amg_alpha_theta);
1067:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS: vector Poisson solver\n");
1068:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG coarsening type %d\n",jac->as_amg_beta_opts[0]);
1069:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG levels of aggressive coarsening %d\n",jac->as_amg_beta_opts[1]);
1070:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG relaxation type %d\n",jac->as_amg_beta_opts[2]);
1071:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG interpolation type %d\n",jac->as_amg_beta_opts[3]);
1072:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG max nonzero elements in interpolation rows %d\n",jac->as_amg_beta_opts[4]);
1073:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS:     boomerAMG strength threshold %g\n",jac->as_amg_beta_theta);
1074:   }
1075:   return(0);
1076: }

1078: static PetscErrorCode PCHYPRESetDiscreteGradient_HYPRE(PC pc, Mat G)
1079: {
1080:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1081:   PetscBool      ishypre;

1085:   PetscObjectTypeCompare((PetscObject)G,MATHYPRE,&ishypre);
1086:   if (ishypre) {
1087:     PetscObjectReference((PetscObject)G);
1088:     MatDestroy(&jac->G);
1089:     jac->G = G;
1090:   } else {
1091:     MatReuse reuse = MAT_INITIAL_MATRIX;
1092:     if (jac->G) reuse = MAT_REUSE_MATRIX;
1093:     MatConvert(G,MATHYPRE,reuse,&jac->G);
1094:   }
1095:   return(0);
1096: }

1098: /*@
1099:  PCHYPRESetDiscreteGradient - Set discrete gradient matrix

1101:    Collective on PC

1103:    Input Parameters:
1104: +  pc - the preconditioning context
1105: -  G - the discrete gradient

1107:    Level: intermediate

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

1112: .seealso:
1113: @*/
1114: PetscErrorCode PCHYPRESetDiscreteGradient(PC pc, Mat G)
1115: {

1122:   PetscTryMethod(pc,"PCHYPRESetDiscreteGradient_C",(PC,Mat),(pc,G));
1123:   return(0);
1124: }

1126: static PetscErrorCode PCHYPRESetDiscreteCurl_HYPRE(PC pc, Mat C)
1127: {
1128:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1129:   PetscBool      ishypre;

1133:   PetscObjectTypeCompare((PetscObject)C,MATHYPRE,&ishypre);
1134:   if (ishypre) {
1135:     PetscObjectReference((PetscObject)C);
1136:     MatDestroy(&jac->C);
1137:     jac->C = C;
1138:   } else {
1139:     MatReuse reuse = MAT_INITIAL_MATRIX;
1140:     if (jac->C) reuse = MAT_REUSE_MATRIX;
1141:     MatConvert(C,MATHYPRE,reuse,&jac->C);
1142:   }
1143:   return(0);
1144: }

1146: /*@
1147:  PCHYPRESetDiscreteCurl - Set discrete curl matrix

1149:    Collective on PC

1151:    Input Parameters:
1152: +  pc - the preconditioning context
1153: -  C - the discrete curl

1155:    Level: intermediate

1157:    Notes: C should have as many rows as the number of faces and as many columns as the number of edges in the mesh
1158:           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

1160: .seealso:
1161: @*/
1162: PetscErrorCode PCHYPRESetDiscreteCurl(PC pc, Mat C)
1163: {

1170:   PetscTryMethod(pc,"PCHYPRESetDiscreteCurl_C",(PC,Mat),(pc,C));
1171:   return(0);
1172: }

1174: static PetscErrorCode PCHYPRESetPoissonMatrix_HYPRE(PC pc, Mat A, PetscBool isalpha)
1175: {
1176:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1177:   PetscBool      ishypre;

1181:   PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);
1182:   if (ishypre) {
1183:     if (isalpha) {
1184:       PetscObjectReference((PetscObject)A);
1185:       MatDestroy(&jac->alpha_Poisson);
1186:       jac->alpha_Poisson = A;
1187:     } else {
1188:       if (A) {
1189:         PetscObjectReference((PetscObject)A);
1190:       } else {
1191:         jac->ams_beta_is_zero = PETSC_TRUE;
1192:       }
1193:       MatDestroy(&jac->beta_Poisson);
1194:       jac->beta_Poisson = A;
1195:     }
1196:   } else {
1197:     if (isalpha) {
1198:       MatReuse reuse = MAT_INITIAL_MATRIX;
1199:       if (jac->alpha_Poisson) reuse = MAT_REUSE_MATRIX;
1200:       MatConvert(A,MATHYPRE,reuse,&jac->alpha_Poisson);
1201:     } else {
1202:       if (A) {
1203:         MatReuse reuse = MAT_INITIAL_MATRIX;
1204:         if (jac->beta_Poisson) reuse = MAT_REUSE_MATRIX;
1205:         MatConvert(A,MATHYPRE,reuse,&jac->beta_Poisson);
1206:       } else {
1207:         MatDestroy(&jac->beta_Poisson);
1208:         jac->ams_beta_is_zero = PETSC_TRUE;
1209:       }
1210:     }
1211:   }
1212:   return(0);
1213: }

1215: /*@
1216:  PCHYPRESetAlphaPoissonMatrix - Set vector Poisson matrix

1218:    Collective on PC

1220:    Input Parameters:
1221: +  pc - the preconditioning context
1222: -  A - the matrix

1224:    Level: intermediate

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

1228: .seealso:
1229: @*/
1230: PetscErrorCode PCHYPRESetAlphaPoissonMatrix(PC pc, Mat A)
1231: {

1238:   PetscTryMethod(pc,"PCHYPRESetPoissonMatrix_C",(PC,Mat,PetscBool),(pc,A,PETSC_TRUE));
1239:   return(0);
1240: }

1242: /*@
1243:  PCHYPRESetBetaPoissonMatrix - Set Poisson matrix

1245:    Collective on PC

1247:    Input Parameters:
1248: +  pc - the preconditioning context
1249: -  A - the matrix

1251:    Level: intermediate

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

1256: .seealso:
1257: @*/
1258: PetscErrorCode PCHYPRESetBetaPoissonMatrix(PC pc, Mat A)
1259: {

1264:   if (A) {
1267:   }
1268:   PetscTryMethod(pc,"PCHYPRESetPoissonMatrix_C",(PC,Mat,PetscBool),(pc,A,PETSC_FALSE));
1269:   return(0);
1270: }

1272: static PetscErrorCode PCHYPRESetEdgeConstantVectors_HYPRE(PC pc,Vec ozz, Vec zoz, Vec zzo)
1273: {
1274:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
1275:   PetscErrorCode     ierr;

1278:   /* throw away any vector if already set */
1279:   if (jac->constants[0]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[0]));
1280:   if (jac->constants[1]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[1]));
1281:   if (jac->constants[2]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[2]));
1282:   jac->constants[0] = NULL;
1283:   jac->constants[1] = NULL;
1284:   jac->constants[2] = NULL;
1285:   VecHYPRE_IJVectorCreate(ozz,&jac->constants[0]);
1286:   VecHYPRE_IJVectorCopy(ozz,jac->constants[0]);
1287:   VecHYPRE_IJVectorCreate(zoz,&jac->constants[1]);
1288:   VecHYPRE_IJVectorCopy(zoz,jac->constants[1]);
1289:   jac->dim = 2;
1290:   if (zzo) {
1291:     VecHYPRE_IJVectorCreate(zzo,&jac->constants[2]);
1292:     VecHYPRE_IJVectorCopy(zzo,jac->constants[2]);
1293:     jac->dim++;
1294:   }
1295:   return(0);
1296: }

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

1301:    Collective on PC

1303:    Input Parameters:
1304: +  pc - the preconditioning context
1305: -  ozz - vector representing (1,0,0) (or (1,0) in 2D)
1306: -  zoz - vector representing (0,1,0) (or (0,1) in 2D)
1307: -  zzo - vector representing (0,0,1) (use NULL in 2D)

1309:    Level: intermediate

1311:    Notes:

1313: .seealso:
1314: @*/
1315: PetscErrorCode PCHYPRESetEdgeConstantVectors(PC pc, Vec ozz, Vec zoz, Vec zzo)
1316: {

1327:   PetscTryMethod(pc,"PCHYPRESetEdgeConstantVectors_C",(PC,Vec,Vec,Vec),(pc,ozz,zoz,zzo));
1328:   return(0);
1329: }

1331: static PetscErrorCode PCSetCoordinates_HYPRE(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords)
1332: {
1333:   PC_HYPRE        *jac = (PC_HYPRE*)pc->data;
1334:   Vec             tv;
1335:   PetscInt        i;
1336:   PetscErrorCode  ierr;

1339:   /* throw away any coordinate vector if already set */
1340:   if (jac->coords[0]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[0]));
1341:   if (jac->coords[1]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[1]));
1342:   if (jac->coords[2]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[2]));
1343:   jac->dim = dim;

1345:   /* compute IJ vector for coordinates */
1346:   VecCreate(PetscObjectComm((PetscObject)pc),&tv);
1347:   VecSetType(tv,VECSTANDARD);
1348:   VecSetSizes(tv,nloc,PETSC_DECIDE);
1349:   for (i=0;i<dim;i++) {
1350:     PetscScalar *array;
1351:     PetscInt    j;

1353:     VecHYPRE_IJVectorCreate(tv,&jac->coords[i]);
1354:     VecGetArray(tv,&array);
1355:     for (j=0;j<nloc;j++) {
1356:       array[j] = coords[j*dim+i];
1357:     }
1358:     PetscStackCallStandard(HYPRE_IJVectorSetValues,(jac->coords[i],nloc,NULL,array));
1359:     PetscStackCallStandard(HYPRE_IJVectorAssemble,(jac->coords[i]));
1360:     VecRestoreArray(tv,&array);
1361:   }
1362:   VecDestroy(&tv);
1363:   return(0);
1364: }

1366: /* ---------------------------------------------------------------------------------*/

1368: static PetscErrorCode  PCHYPREGetType_HYPRE(PC pc,const char *name[])
1369: {
1370:   PC_HYPRE *jac = (PC_HYPRE*)pc->data;

1373:   *name = jac->hypre_type;
1374:   return(0);
1375: }

1377: static PetscErrorCode  PCHYPRESetType_HYPRE(PC pc,const char name[])
1378: {
1379:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1381:   PetscBool      flag;

1384:   if (jac->hypre_type) {
1385:     PetscStrcmp(jac->hypre_type,name,&flag);
1386:     if (!flag) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
1387:     return(0);
1388:   } else {
1389:     PetscStrallocpy(name, &jac->hypre_type);
1390:   }

1392:   jac->maxiter         = PETSC_DEFAULT;
1393:   jac->tol             = PETSC_DEFAULT;
1394:   jac->printstatistics = PetscLogPrintInfo;

1396:   PetscStrcmp("pilut",jac->hypre_type,&flag);
1397:   if (flag) {
1398:     PetscStackCallStandard(HYPRE_ParCSRPilutCreate,(jac->comm_hypre,&jac->hsolver));
1399:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
1400:     pc->ops->view           = PCView_HYPRE_Pilut;
1401:     jac->destroy            = HYPRE_ParCSRPilutDestroy;
1402:     jac->setup              = HYPRE_ParCSRPilutSetup;
1403:     jac->solve              = HYPRE_ParCSRPilutSolve;
1404:     jac->factorrowsize      = PETSC_DEFAULT;
1405:     return(0);
1406:   }
1407:   PetscStrcmp("parasails",jac->hypre_type,&flag);
1408:   if (flag) {
1409:     PetscStackCallStandard(HYPRE_ParaSailsCreate,(jac->comm_hypre,&jac->hsolver));
1410:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
1411:     pc->ops->view           = PCView_HYPRE_ParaSails;
1412:     jac->destroy            = HYPRE_ParaSailsDestroy;
1413:     jac->setup              = HYPRE_ParaSailsSetup;
1414:     jac->solve              = HYPRE_ParaSailsSolve;
1415:     /* initialize */
1416:     jac->nlevels    = 1;
1417:     jac->threshhold = .1;
1418:     jac->filter     = .1;
1419:     jac->loadbal    = 0;
1420:     if (PetscLogPrintInfo) jac->logging = (int) PETSC_TRUE;
1421:     else jac->logging = (int) PETSC_FALSE;

1423:     jac->ruse = (int) PETSC_FALSE;
1424:     jac->symt = 0;
1425:     PetscStackCallStandard(HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels));
1426:     PetscStackCallStandard(HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter));
1427:     PetscStackCallStandard(HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal));
1428:     PetscStackCallStandard(HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging));
1429:     PetscStackCallStandard(HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse));
1430:     PetscStackCallStandard(HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
1431:     return(0);
1432:   }
1433:   PetscStrcmp("boomeramg",jac->hypre_type,&flag);
1434:   if (flag) {
1435:     HYPRE_BoomerAMGCreate(&jac->hsolver);
1436:     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
1437:     pc->ops->view            = PCView_HYPRE_BoomerAMG;
1438:     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
1439:     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
1440:     jac->destroy             = HYPRE_BoomerAMGDestroy;
1441:     jac->setup               = HYPRE_BoomerAMGSetup;
1442:     jac->solve               = HYPRE_BoomerAMGSolve;
1443:     jac->applyrichardson     = PETSC_FALSE;
1444:     /* these defaults match the hypre defaults */
1445:     jac->cycletype        = 1;
1446:     jac->maxlevels        = 25;
1447:     jac->maxiter          = 1;
1448:     jac->tol              = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */
1449:     jac->truncfactor      = 0.0;
1450:     jac->strongthreshold  = .25;
1451:     jac->maxrowsum        = .9;
1452:     jac->coarsentype      = 6;
1453:     jac->measuretype      = 0;
1454:     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = 1;
1455:     jac->smoothtype       = -1; /* Not set by default */
1456:     jac->smoothnumlevels  = 25;
1457:     jac->eu_level         = 0;
1458:     jac->eu_droptolerance = 0;
1459:     jac->eu_bj            = 0;
1460:     jac->relaxtype[0]     = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */
1461:     jac->relaxtype[2]     = 9; /*G.E. */
1462:     jac->relaxweight      = 1.0;
1463:     jac->outerrelaxweight = 1.0;
1464:     jac->relaxorder       = 1;
1465:     jac->interptype       = 0;
1466:     jac->agg_nl           = 0;
1467:     jac->pmax             = 0;
1468:     jac->truncfactor      = 0.0;
1469:     jac->agg_num_paths    = 1;

1471:     jac->nodal_coarsen      = 0;
1472:     jac->nodal_relax        = PETSC_FALSE;
1473:     jac->nodal_relax_levels = 1;
1474:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
1475:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
1476:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
1477:     PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
1478:     PetscStackCallStandard(HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
1479:     PetscStackCallStandard(HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
1480:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
1481:     PetscStackCallStandard(HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
1482:     PetscStackCallStandard(HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
1483:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
1484:     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
1485:     PetscStackCallStandard(HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
1486:     PetscStackCallStandard(HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
1487:     PetscStackCallStandard(HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
1488:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, jac->relaxtype[0]));  /*defaults coarse to 9*/
1489:     PetscStackCallStandard(HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver, jac->gridsweeps[0])); /*defaults coarse to 1 */
1490:     return(0);
1491:   }
1492:   PetscStrcmp("ams",jac->hypre_type,&flag);
1493:   if (flag) {
1494:     HYPRE_AMSCreate(&jac->hsolver);
1495:     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_AMS;
1496:     pc->ops->view            = PCView_HYPRE_AMS;
1497:     jac->destroy             = HYPRE_AMSDestroy;
1498:     jac->setup               = HYPRE_AMSSetup;
1499:     jac->solve               = HYPRE_AMSSolve;
1500:     jac->coords[0]           = NULL;
1501:     jac->coords[1]           = NULL;
1502:     jac->coords[2]           = NULL;
1503:     /* solver parameters: these are borrowed from mfem package, and they are not the default values from HYPRE AMS */
1504:     jac->as_print           = 0;
1505:     jac->as_max_iter        = 1; /* used as a preconditioner */
1506:     jac->as_tol             = 0.; /* used as a preconditioner */
1507:     jac->ams_cycle_type     = 13;
1508:     /* Smoothing options */
1509:     jac->as_relax_type      = 2;
1510:     jac->as_relax_times     = 1;
1511:     jac->as_relax_weight    = 1.0;
1512:     jac->as_omega           = 1.0;
1513:     /* Vector valued Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
1514:     jac->as_amg_alpha_opts[0] = 10;
1515:     jac->as_amg_alpha_opts[1] = 1;
1516:     jac->as_amg_alpha_opts[2] = 6;
1517:     jac->as_amg_alpha_opts[3] = 6;
1518:     jac->as_amg_alpha_opts[4] = 4;
1519:     jac->as_amg_alpha_theta   = 0.25;
1520:     /* Scalar Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
1521:     jac->as_amg_beta_opts[0] = 10;
1522:     jac->as_amg_beta_opts[1] = 1;
1523:     jac->as_amg_beta_opts[2] = 6;
1524:     jac->as_amg_beta_opts[3] = 6;
1525:     jac->as_amg_beta_opts[4] = 4;
1526:     jac->as_amg_beta_theta   = 0.25;
1527:     PetscStackCallStandard(HYPRE_AMSSetPrintLevel,(jac->hsolver,jac->as_print));
1528:     PetscStackCallStandard(HYPRE_AMSSetMaxIter,(jac->hsolver,jac->as_max_iter));
1529:     PetscStackCallStandard(HYPRE_AMSSetCycleType,(jac->hsolver,jac->ams_cycle_type));
1530:     PetscStackCallStandard(HYPRE_AMSSetTol,(jac->hsolver,jac->as_tol));
1531:     PetscStackCallStandard(HYPRE_AMSSetSmoothingOptions,(jac->hsolver,jac->as_relax_type,
1532:                                                                       jac->as_relax_times,
1533:                                                                       jac->as_relax_weight,
1534:                                                                       jac->as_omega));
1535:     PetscStackCallStandard(HYPRE_AMSSetAlphaAMGOptions,(jac->hsolver,jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
1536:                                                                      jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
1537:                                                                      jac->as_amg_alpha_opts[2],       /* AMG relax_type */
1538:                                                                      jac->as_amg_alpha_theta,
1539:                                                                      jac->as_amg_alpha_opts[3],       /* AMG interp_type */
1540:                                                                      jac->as_amg_alpha_opts[4]));     /* AMG Pmax */
1541:     PetscStackCallStandard(HYPRE_AMSSetBetaAMGOptions,(jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
1542:                                                                     jac->as_amg_beta_opts[1],       /* AMG agg_levels */
1543:                                                                     jac->as_amg_beta_opts[2],       /* AMG relax_type */
1544:                                                                     jac->as_amg_beta_theta,
1545:                                                                     jac->as_amg_beta_opts[3],       /* AMG interp_type */
1546:                                                                     jac->as_amg_beta_opts[4]));     /* AMG Pmax */
1547:     /* Zero conductivity */
1548:     jac->ams_beta_is_zero      = PETSC_FALSE;
1549:     jac->ams_beta_is_zero_part = PETSC_FALSE;
1550:     return(0);
1551:   }
1552:   PetscStrcmp("ads",jac->hypre_type,&flag);
1553:   if (flag) {
1554:     HYPRE_ADSCreate(&jac->hsolver);
1555:     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_ADS;
1556:     pc->ops->view            = PCView_HYPRE_ADS;
1557:     jac->destroy             = HYPRE_ADSDestroy;
1558:     jac->setup               = HYPRE_ADSSetup;
1559:     jac->solve               = HYPRE_ADSSolve;
1560:     jac->coords[0]           = NULL;
1561:     jac->coords[1]           = NULL;
1562:     jac->coords[2]           = NULL;
1563:     /* solver parameters: these are borrowed from mfem package, and they are not the default values from HYPRE ADS */
1564:     jac->as_print           = 0;
1565:     jac->as_max_iter        = 1; /* used as a preconditioner */
1566:     jac->as_tol             = 0.; /* used as a preconditioner */
1567:     jac->ads_cycle_type     = 13;
1568:     /* Smoothing options */
1569:     jac->as_relax_type      = 2;
1570:     jac->as_relax_times     = 1;
1571:     jac->as_relax_weight    = 1.0;
1572:     jac->as_omega           = 1.0;
1573:     /* AMS solver parameters: cycle_type, coarsen type, agg_levels, relax_type, interp_type, Pmax */
1574:     jac->ams_cycle_type       = 14;
1575:     jac->as_amg_alpha_opts[0] = 10;
1576:     jac->as_amg_alpha_opts[1] = 1;
1577:     jac->as_amg_alpha_opts[2] = 6;
1578:     jac->as_amg_alpha_opts[3] = 6;
1579:     jac->as_amg_alpha_opts[4] = 4;
1580:     jac->as_amg_alpha_theta   = 0.25;
1581:     /* Vector Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
1582:     jac->as_amg_beta_opts[0] = 10;
1583:     jac->as_amg_beta_opts[1] = 1;
1584:     jac->as_amg_beta_opts[2] = 6;
1585:     jac->as_amg_beta_opts[3] = 6;
1586:     jac->as_amg_beta_opts[4] = 4;
1587:     jac->as_amg_beta_theta   = 0.25;
1588:     PetscStackCallStandard(HYPRE_ADSSetPrintLevel,(jac->hsolver,jac->as_print));
1589:     PetscStackCallStandard(HYPRE_ADSSetMaxIter,(jac->hsolver,jac->as_max_iter));
1590:     PetscStackCallStandard(HYPRE_ADSSetCycleType,(jac->hsolver,jac->ams_cycle_type));
1591:     PetscStackCallStandard(HYPRE_ADSSetTol,(jac->hsolver,jac->as_tol));
1592:     PetscStackCallStandard(HYPRE_ADSSetSmoothingOptions,(jac->hsolver,jac->as_relax_type,
1593:                                                                       jac->as_relax_times,
1594:                                                                       jac->as_relax_weight,
1595:                                                                       jac->as_omega));
1596:     PetscStackCallStandard(HYPRE_ADSSetAMSOptions,(jac->hsolver,jac->ams_cycle_type,             /* AMG coarsen type */
1597:                                                                 jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
1598:                                                                 jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
1599:                                                                 jac->as_amg_alpha_opts[2],       /* AMG relax_type */
1600:                                                                 jac->as_amg_alpha_theta,
1601:                                                                 jac->as_amg_alpha_opts[3],       /* AMG interp_type */
1602:                                                                 jac->as_amg_alpha_opts[4]));     /* AMG Pmax */
1603:     PetscStackCallStandard(HYPRE_ADSSetAMGOptions,(jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
1604:                                                                 jac->as_amg_beta_opts[1],       /* AMG agg_levels */
1605:                                                                 jac->as_amg_beta_opts[2],       /* AMG relax_type */
1606:                                                                 jac->as_amg_beta_theta,
1607:                                                                 jac->as_amg_beta_opts[3],       /* AMG interp_type */
1608:                                                                 jac->as_amg_beta_opts[4]));     /* AMG Pmax */
1609:     return(0);
1610:   }
1611:   PetscFree(jac->hypre_type);

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

1618: /*
1619:     It only gets here if the HYPRE type has not been set before the call to
1620:    ...SetFromOptions() which actually is most of the time
1621: */
1622: static PetscErrorCode PCSetFromOptions_HYPRE(PetscOptionItems *PetscOptionsObject,PC pc)
1623: {
1625:   PetscInt       indx;
1626:   const char     *type[] = {"pilut","parasails","boomeramg","ams","ads"};
1627:   PetscBool      flg;

1630:   PetscOptionsHead(PetscOptionsObject,"HYPRE preconditioner options");
1631:   PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"boomeramg",&indx,&flg);
1632:   if (flg) {
1633:     PCHYPRESetType_HYPRE(pc,type[indx]);
1634:   } else {
1635:     PCHYPRESetType_HYPRE(pc,"boomeramg");
1636:   }
1637:   if (pc->ops->setfromoptions) {
1638:     pc->ops->setfromoptions(PetscOptionsObject,pc);
1639:   }
1640:   PetscOptionsTail();
1641:   return(0);
1642: }

1644: /*@C
1645:      PCHYPRESetType - Sets which hypre preconditioner you wish to use

1647:    Input Parameters:
1648: +     pc - the preconditioner context
1649: -     name - either  pilut, parasails, boomeramg, ams, ads

1651:    Options Database Keys:
1652:    -pc_hypre_type - One of pilut, parasails, boomeramg, ams, ads

1654:    Level: intermediate

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

1659: @*/
1660: PetscErrorCode  PCHYPRESetType(PC pc,const char name[])
1661: {

1667:   PetscTryMethod(pc,"PCHYPRESetType_C",(PC,const char[]),(pc,name));
1668:   return(0);
1669: }

1671: /*@C
1672:      PCHYPREGetType - Gets which hypre preconditioner you are using

1674:    Input Parameter:
1675: .     pc - the preconditioner context

1677:    Output Parameter:
1678: .     name - either  pilut, parasails, boomeramg, ams, ads

1680:    Level: intermediate

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

1685: @*/
1686: PetscErrorCode  PCHYPREGetType(PC pc,const char *name[])
1687: {

1693:   PetscTryMethod(pc,"PCHYPREGetType_C",(PC,const char*[]),(pc,name));
1694:   return(0);
1695: }

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

1700:    Options Database Keys:
1701: +   -pc_hypre_type - One of pilut, parasails, boomeramg, ams, ads
1702: -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
1703:           preconditioner

1705:    Level: intermediate

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

1711:           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_tol refer to the number of iterations
1712:           (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if
1713:           -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner
1714:           (-pc_hypre_boomeramg_tol should be set to 0.0 - the default - to strictly use a fixed number of
1715:           iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations
1716:           and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10
1717:           then AT MOST twenty V-cycles of boomeramg will be called.

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

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

1729:           MatSetNearNullSpace() - if you provide a near null space to your matrix it is ignored by hypre UNLESS you also use
1730:           the two options:
1731: +   -pc_hypre_boomeramg_nodal_coarsen <n> - where n is from 1 to 6 (see HYPRE_BOOMERAMGSetNodal())
1732: -   -pc_hypre_boomeramg_vec_interp_variant <v> where v is from 1 to 3 (see HYPRE_BoomerAMGSetInterpVecVariant())

1734:           Depending on the linear system you may see the same or different convergence depending on the values you use.

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

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

1741: M*/

1743: PETSC_EXTERN PetscErrorCode PCCreate_HYPRE(PC pc)
1744: {
1745:   PC_HYPRE       *jac;

1749:   PetscNewLog(pc,&jac);

1751:   pc->data                = jac;
1752:   pc->ops->reset          = PCReset_HYPRE;
1753:   pc->ops->destroy        = PCDestroy_HYPRE;
1754:   pc->ops->setfromoptions = PCSetFromOptions_HYPRE;
1755:   pc->ops->setup          = PCSetUp_HYPRE;
1756:   pc->ops->apply          = PCApply_HYPRE;
1757:   jac->comm_hypre         = MPI_COMM_NULL;
1758:   jac->coords[0]          = NULL;
1759:   jac->coords[1]          = NULL;
1760:   jac->coords[2]          = NULL;
1761:   jac->constants[0]       = NULL;
1762:   jac->constants[1]       = NULL;
1763:   jac->constants[2]       = NULL;
1764:   jac->hmnull             = NULL;
1765:   jac->n_hmnull           = 0;
1766:   /* duplicate communicator for hypre */
1767:   MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(jac->comm_hypre));
1768:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetType_C",PCHYPRESetType_HYPRE);
1769:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPREGetType_C",PCHYPREGetType_HYPRE);
1770:   PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_HYPRE);
1771:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteGradient_C",PCHYPRESetDiscreteGradient_HYPRE);
1772:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteCurl_C",PCHYPRESetDiscreteCurl_HYPRE);
1773:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetEdgeConstantVectors_C",PCHYPRESetEdgeConstantVectors_HYPRE);
1774:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetPoissonMatrix_C",PCHYPRESetPoissonMatrix_HYPRE);
1775:   return(0);
1776: }

1778: /* ---------------------------------------------------------------------------------------------------------------------------------*/

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

1784:   /* keep copy of PFMG options used so may view them */
1785:   PetscInt its;
1786:   double   tol;
1787:   PetscInt relax_type;
1788:   PetscInt rap_type;
1789:   PetscInt num_pre_relax,num_post_relax;
1790:   PetscInt max_levels;
1791: } PC_PFMG;

1793: PetscErrorCode PCDestroy_PFMG(PC pc)
1794: {
1796:   PC_PFMG        *ex = (PC_PFMG*) pc->data;

1799:   if (ex->hsolver) PetscStackCallStandard(HYPRE_StructPFMGDestroy,(ex->hsolver));
1800:   MPI_Comm_free(&ex->hcomm);
1801:   PetscFree(pc->data);
1802:   return(0);
1803: }

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

1808: PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer)
1809: {
1811:   PetscBool      iascii;
1812:   PC_PFMG        *ex = (PC_PFMG*) pc->data;

1815:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1816:   if (iascii) {
1817:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG preconditioning\n");
1818:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max iterations %d\n",ex->its);
1819:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: tolerance %g\n",ex->tol);
1820:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);
1821:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: RAP type %s\n",PFMGRAPType[ex->rap_type]);
1822:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);
1823:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max levels %d\n",ex->max_levels);
1824:   }
1825:   return(0);
1826: }


1829: PetscErrorCode PCSetFromOptions_PFMG(PetscOptionItems *PetscOptionsObject,PC pc)
1830: {
1832:   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1833:   PetscBool      flg = PETSC_FALSE;

1836:   PetscOptionsHead(PetscOptionsObject,"PFMG options");
1837:   PetscOptionsBool("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",flg,&flg,NULL);
1838:   if (flg) {
1839:     int level=3;
1840:     PetscStackCallStandard(HYPRE_StructPFMGSetPrintLevel,(ex->hsolver,level));
1841:   }
1842:   PetscOptionsInt("-pc_pfmg_its","Number of iterations of PFMG to use as preconditioner","HYPRE_StructPFMGSetMaxIter",ex->its,&ex->its,NULL);
1843:   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(ex->hsolver,ex->its));
1844:   PetscOptionsInt("-pc_pfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_StructPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,NULL);
1845:   PetscStackCallStandard(HYPRE_StructPFMGSetNumPreRelax,(ex->hsolver,ex->num_pre_relax));
1846:   PetscOptionsInt("-pc_pfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_StructPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,NULL);
1847:   PetscStackCallStandard(HYPRE_StructPFMGSetNumPostRelax,(ex->hsolver,ex->num_post_relax));

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

1852:   PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,NULL);
1853:   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(ex->hsolver,ex->tol));
1854:   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);
1855:   PetscStackCallStandard(HYPRE_StructPFMGSetRelaxType,(ex->hsolver, ex->relax_type));
1856:   PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,ALEN(PFMGRAPType),PFMGRAPType[ex->rap_type],&ex->rap_type,NULL);
1857:   PetscStackCallStandard(HYPRE_StructPFMGSetRAPType,(ex->hsolver, ex->rap_type));
1858:   PetscOptionsTail();
1859:   return(0);
1860: }

1862: PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y)
1863: {
1864:   PetscErrorCode    ierr;
1865:   PC_PFMG           *ex = (PC_PFMG*) pc->data;
1866:   PetscScalar       *yy;
1867:   const PetscScalar *xx;
1868:   PetscInt          ilower[3],iupper[3];
1869:   Mat_HYPREStruct   *mx = (Mat_HYPREStruct*)(pc->pmat->data);

1872:   PetscCitationsRegister(hypreCitation,&cite);
1873:   DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);
1874:   iupper[0] += ilower[0] - 1;
1875:   iupper[1] += ilower[1] - 1;
1876:   iupper[2] += ilower[2] - 1;

1878:   /* copy x values over to hypre */
1879:   PetscStackCallStandard(HYPRE_StructVectorSetConstantValues,(mx->hb,0.0));
1880:   VecGetArrayRead(x,&xx);
1881:   PetscStackCallStandard(HYPRE_StructVectorSetBoxValues,(mx->hb,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,(PetscScalar*)xx));
1882:   VecRestoreArrayRead(x,&xx);
1883:   PetscStackCallStandard(HYPRE_StructVectorAssemble,(mx->hb));
1884:   PetscStackCallStandard(HYPRE_StructPFMGSolve,(ex->hsolver,mx->hmat,mx->hb,mx->hx));

1886:   /* copy solution values back to PETSc */
1887:   VecGetArray(y,&yy);
1888:   PetscStackCallStandard(HYPRE_StructVectorGetBoxValues,(mx->hx,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,yy));
1889:   VecRestoreArray(y,&yy);
1890:   return(0);
1891: }

1893: 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)
1894: {
1895:   PC_PFMG        *jac = (PC_PFMG*)pc->data;
1897:   PetscInt       oits;

1900:   PetscCitationsRegister(hypreCitation,&cite);
1901:   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(jac->hsolver,its*jac->its));
1902:   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(jac->hsolver,rtol));

1904:   PCApply_PFMG(pc,b,y);
1905:   PetscStackCallStandard(HYPRE_StructPFMGGetNumIterations,(jac->hsolver,(HYPRE_Int *)&oits));
1906:   *outits = oits;
1907:   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
1908:   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
1909:   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(jac->hsolver,jac->tol));
1910:   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(jac->hsolver,jac->its));
1911:   return(0);
1912: }


1915: PetscErrorCode PCSetUp_PFMG(PC pc)
1916: {
1917:   PetscErrorCode  ierr;
1918:   PC_PFMG         *ex = (PC_PFMG*) pc->data;
1919:   Mat_HYPREStruct *mx = (Mat_HYPREStruct*)(pc->pmat->data);
1920:   PetscBool       flg;

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

1926:   /* create the hypre solver object and set its information */
1927:   if (ex->hsolver) PetscStackCallStandard(HYPRE_StructPFMGDestroy,(ex->hsolver));
1928:   PetscStackCallStandard(HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
1929:   PetscStackCallStandard(HYPRE_StructPFMGSetup,(ex->hsolver,mx->hmat,mx->hb,mx->hx));
1930:   PetscStackCallStandard(HYPRE_StructPFMGSetZeroGuess,(ex->hsolver));
1931:   return(0);
1932: }


1935: /*MC
1936:      PCPFMG - the hypre PFMG multigrid solver

1938:    Level: advanced

1940:    Options Database:
1941: + -pc_pfmg_its <its> number of iterations of PFMG to use as preconditioner
1942: . -pc_pfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
1943: . -pc_pfmg_num_post_relax <steps> number of smoothing steps after coarse grid
1944: . -pc_pfmg_tol <tol> tolerance of PFMG
1945: . -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
1946: - -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin

1948:    Notes:  This is for CELL-centered descretizations

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

1953: .seealso:  PCMG, MATHYPRESTRUCT
1954: M*/

1956: PETSC_EXTERN PetscErrorCode PCCreate_PFMG(PC pc)
1957: {
1959:   PC_PFMG        *ex;

1962:   PetscNew(&ex); \
1963:   pc->data = ex;

1965:   ex->its            = 1;
1966:   ex->tol            = 1.e-8;
1967:   ex->relax_type     = 1;
1968:   ex->rap_type       = 0;
1969:   ex->num_pre_relax  = 1;
1970:   ex->num_post_relax = 1;
1971:   ex->max_levels     = 0;

1973:   pc->ops->setfromoptions  = PCSetFromOptions_PFMG;
1974:   pc->ops->view            = PCView_PFMG;
1975:   pc->ops->destroy         = PCDestroy_PFMG;
1976:   pc->ops->apply           = PCApply_PFMG;
1977:   pc->ops->applyrichardson = PCApplyRichardson_PFMG;
1978:   pc->ops->setup           = PCSetUp_PFMG;

1980:   MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ex->hcomm));
1981:   PetscStackCallStandard(HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
1982:   return(0);
1983: }

1985: /* ---------------------------------------------------------------------------------------------------------------------------------------------------*/

1987: /* we know we are working with a HYPRE_SStructMatrix */
1988: typedef struct {
1989:   MPI_Comm            hcomm;       /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */
1990:   HYPRE_SStructSolver ss_solver;

1992:   /* keep copy of SYSPFMG options used so may view them */
1993:   PetscInt its;
1994:   double   tol;
1995:   PetscInt relax_type;
1996:   PetscInt num_pre_relax,num_post_relax;
1997: } PC_SysPFMG;

1999: PetscErrorCode PCDestroy_SysPFMG(PC pc)
2000: {
2002:   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;

2005:   if (ex->ss_solver) PetscStackCallStandard(HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));
2006:   MPI_Comm_free(&ex->hcomm);
2007:   PetscFree(pc->data);
2008:   return(0);
2009: }

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

2013: PetscErrorCode PCView_SysPFMG(PC pc,PetscViewer viewer)
2014: {
2016:   PetscBool      iascii;
2017:   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;

2020:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
2021:   if (iascii) {
2022:     PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG preconditioning\n");
2023:     PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: max iterations %d\n",ex->its);
2024:     PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: tolerance %g\n",ex->tol);
2025:     PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);
2026:     PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);
2027:   }
2028:   return(0);
2029: }


2032: PetscErrorCode PCSetFromOptions_SysPFMG(PetscOptionItems *PetscOptionsObject,PC pc)
2033: {
2035:   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;
2036:   PetscBool      flg = PETSC_FALSE;

2039:   PetscOptionsHead(PetscOptionsObject,"SysPFMG options");
2040:   PetscOptionsBool("-pc_syspfmg_print_statistics","Print statistics","HYPRE_SStructSysPFMGSetPrintLevel",flg,&flg,NULL);
2041:   if (flg) {
2042:     int level=3;
2043:     PetscStackCallStandard(HYPRE_SStructSysPFMGSetPrintLevel,(ex->ss_solver,level));
2044:   }
2045:   PetscOptionsInt("-pc_syspfmg_its","Number of iterations of SysPFMG to use as preconditioner","HYPRE_SStructSysPFMGSetMaxIter",ex->its,&ex->its,NULL);
2046:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(ex->ss_solver,ex->its));
2047:   PetscOptionsInt("-pc_syspfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_SStructSysPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,NULL);
2048:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetNumPreRelax,(ex->ss_solver,ex->num_pre_relax));
2049:   PetscOptionsInt("-pc_syspfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_SStructSysPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,NULL);
2050:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetNumPostRelax,(ex->ss_solver,ex->num_post_relax));

2052:   PetscOptionsReal("-pc_syspfmg_tol","Tolerance of SysPFMG","HYPRE_SStructSysPFMGSetTol",ex->tol,&ex->tol,NULL);
2053:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(ex->ss_solver,ex->tol));
2054:   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);
2055:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetRelaxType,(ex->ss_solver, ex->relax_type));
2056:   PetscOptionsTail();
2057:   return(0);
2058: }

2060: PetscErrorCode PCApply_SysPFMG(PC pc,Vec x,Vec y)
2061: {
2062:   PetscErrorCode    ierr;
2063:   PC_SysPFMG        *ex = (PC_SysPFMG*) pc->data;
2064:   PetscScalar       *yy;
2065:   const PetscScalar *xx;
2066:   PetscInt          ilower[3],iupper[3];
2067:   Mat_HYPRESStruct  *mx     = (Mat_HYPRESStruct*)(pc->pmat->data);
2068:   PetscInt          ordering= mx->dofs_order;
2069:   PetscInt          nvars   = mx->nvars;
2070:   PetscInt          part    = 0;
2071:   PetscInt          size;
2072:   PetscInt          i;

2075:   PetscCitationsRegister(hypreCitation,&cite);
2076:   DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);
2077:   iupper[0] += ilower[0] - 1;
2078:   iupper[1] += ilower[1] - 1;
2079:   iupper[2] += ilower[2] - 1;

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

2084:   /* copy x values over to hypre for variable ordering */
2085:   if (ordering) {
2086:     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
2087:     VecGetArrayRead(x,&xx);
2088:     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,(PetscScalar*)xx+(size*i)));
2089:     VecRestoreArrayRead(x,&xx);
2090:     PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
2091:     PetscStackCallStandard(HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
2092:     PetscStackCallStandard(HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));

2094:     /* copy solution values back to PETSc */
2095:     VecGetArray(y,&yy);
2096:     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,yy+(size*i)));
2097:     VecRestoreArray(y,&yy);
2098:   } else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
2099:     PetscScalar *z;
2100:     PetscInt    j, k;

2102:     PetscMalloc1(nvars*size,&z);
2103:     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
2104:     VecGetArrayRead(x,&xx);

2106:     /* transform nodal to hypre's variable ordering for sys_pfmg */
2107:     for (i= 0; i< size; i++) {
2108:       k= i*nvars;
2109:       for (j= 0; j< nvars; j++) z[j*size+i]= xx[k+j];
2110:     }
2111:     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,z+(size*i)));
2112:     VecRestoreArrayRead(x,&xx);
2113:     PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
2114:     PetscStackCallStandard(HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));

2116:     /* copy solution values back to PETSc */
2117:     VecGetArray(y,&yy);
2118:     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,z+(size*i)));
2119:     /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
2120:     for (i= 0; i< size; i++) {
2121:       k= i*nvars;
2122:       for (j= 0; j< nvars; j++) yy[k+j]= z[j*size+i];
2123:     }
2124:     VecRestoreArray(y,&yy);
2125:     PetscFree(z);
2126:   }
2127:   return(0);
2128: }

2130: 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)
2131: {
2132:   PC_SysPFMG     *jac = (PC_SysPFMG*)pc->data;
2134:   PetscInt       oits;

2137:   PetscCitationsRegister(hypreCitation,&cite);
2138:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,its*jac->its));
2139:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,rtol));
2140:   PCApply_SysPFMG(pc,b,y);
2141:   PetscStackCallStandard(HYPRE_SStructSysPFMGGetNumIterations,(jac->ss_solver,(HYPRE_Int *)&oits));
2142:   *outits = oits;
2143:   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
2144:   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
2145:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,jac->tol));
2146:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,jac->its));
2147:   return(0);
2148: }


2151: PetscErrorCode PCSetUp_SysPFMG(PC pc)
2152: {
2153:   PetscErrorCode   ierr;
2154:   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
2155:   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct*)(pc->pmat->data);
2156:   PetscBool        flg;

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

2162:   /* create the hypre sstruct solver object and set its information */
2163:   if (ex->ss_solver) PetscStackCallStandard(HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));
2164:   PetscStackCallStandard(HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
2165:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetZeroGuess,(ex->ss_solver));
2166:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetup,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
2167:   return(0);
2168: }


2171: /*MC
2172:      PCSysPFMG - the hypre SysPFMG multigrid solver

2174:    Level: advanced

2176:    Options Database:
2177: + -pc_syspfmg_its <its> number of iterations of SysPFMG to use as preconditioner
2178: . -pc_syspfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
2179: . -pc_syspfmg_num_post_relax <steps> number of smoothing steps after coarse grid
2180: . -pc_syspfmg_tol <tol> tolerance of SysPFMG
2181: . -pc_syspfmg_relax_type -relaxation type for the up and down cycles, one of Weighted-Jacobi,Red/Black-Gauss-Seidel

2183:    Notes:  This is for CELL-centered descretizations

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

2189: .seealso:  PCMG, MATHYPRESSTRUCT
2190: M*/

2192: PETSC_EXTERN PetscErrorCode PCCreate_SysPFMG(PC pc)
2193: {
2195:   PC_SysPFMG     *ex;

2198:   PetscNew(&ex); \
2199:   pc->data = ex;

2201:   ex->its            = 1;
2202:   ex->tol            = 1.e-8;
2203:   ex->relax_type     = 1;
2204:   ex->num_pre_relax  = 1;
2205:   ex->num_post_relax = 1;

2207:   pc->ops->setfromoptions  = PCSetFromOptions_SysPFMG;
2208:   pc->ops->view            = PCView_SysPFMG;
2209:   pc->ops->destroy         = PCDestroy_SysPFMG;
2210:   pc->ops->apply           = PCApply_SysPFMG;
2211:   pc->ops->applyrichardson = PCApplyRichardson_SysPFMG;
2212:   pc->ops->setup           = PCSetUp_SysPFMG;

2214:   MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ex->hcomm));
2215:   PetscStackCallStandard(HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
2216:   return(0);
2217: }