Actual source code: hypre.c

petsc-3.11.0 2019-03-29
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:   PetscBool nodal_relax;
 77:   PetscInt  nodal_relax_levels;

 79:   PetscInt  nodal_coarsening;
 80:   PetscInt  nodal_coarsening_diag;
 81:   PetscInt  vec_interp_variant;
 82:   PetscInt  vec_interp_qmax;
 83:   PetscBool vec_interp_smooth;
 84:   PetscInt  interp_refine;

 86:   HYPRE_IJVector  *hmnull;
 87:   HYPRE_ParVector *phmnull;  /* near null space passed to hypre */
 88:   PetscInt        n_hmnull;
 89:   Vec             hmnull_constant;
 90:   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 */

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

107:   /* additional data */
108:   Mat G;             /* MatHYPRE */
109:   Mat C;             /* MatHYPRE */
110:   Mat alpha_Poisson; /* MatHYPRE */
111:   Mat beta_Poisson;  /* MatHYPRE */

113:   /* extra information for AMS */
114:   PetscInt       dim; /* geometrical dimension */
115:   HYPRE_IJVector coords[3];
116:   HYPRE_IJVector constants[3];
117:   Mat            RT_PiFull, RT_Pi[3];
118:   Mat            ND_PiFull, ND_Pi[3];
119:   PetscBool      ams_beta_is_zero;
120:   PetscBool      ams_beta_is_zero_part;
121:   PetscInt       ams_proj_freq;
122: } PC_HYPRE;

124: PetscErrorCode PCHYPREGetSolver(PC pc,HYPRE_Solver *hsolver)
125: {
126:   PC_HYPRE *jac = (PC_HYPRE*)pc->data;

129:   *hsolver = jac->hsolver;
130:   return(0);
131: }

133: static PetscErrorCode PCSetUp_HYPRE(PC pc)
134: {
135:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
136:   Mat_HYPRE          *hjac;
137:   HYPRE_ParCSRMatrix hmat;
138:   HYPRE_ParVector    bv,xv;
139:   PetscBool          ishypre;
140:   PetscErrorCode     ierr;

143:   if (!jac->hypre_type) {
144:     PCHYPRESetType(pc,"boomeramg");
145:   }

147:   PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRE,&ishypre);
148:   if (!ishypre) {
149:     MatDestroy(&jac->hpmat);
150:     MatConvert(pc->pmat,MATHYPRE,MAT_INITIAL_MATRIX,&jac->hpmat);
151:   } else {
152:     PetscObjectReference((PetscObject)pc->pmat);
153:     MatDestroy(&jac->hpmat);
154:     jac->hpmat = pc->pmat;
155:   }
156:   hjac = (Mat_HYPRE*)(jac->hpmat->data);

158:   /* special case for BoomerAMG */
159:   if (jac->setup == HYPRE_BoomerAMGSetup) {
160:     MatNullSpace    mnull;
161:     PetscBool       has_const;
162:     PetscInt        bs,nvec,i;
163:     const Vec       *vecs;
164:     PetscScalar     *petscvecarray;

166:     MatGetBlockSize(pc->pmat,&bs);
167:     if (bs > 1) PetscStackCallStandard(HYPRE_BoomerAMGSetNumFunctions,(jac->hsolver,bs));
168:     MatGetNearNullSpace(pc->mat, &mnull);
169:     if (mnull) {
170:       MatNullSpaceGetVecs(mnull, &has_const, &nvec, &vecs);
171:       PetscMalloc1(nvec+1,&jac->hmnull);
172:       PetscMalloc1(nvec+1,&jac->hmnull_hypre_data_array);
173:       PetscMalloc1(nvec+1,&jac->phmnull);
174:       for (i=0; i<nvec; i++) {
175:         VecHYPRE_IJVectorCreate(vecs[i],&jac->hmnull[i]);
176:         VecGetArrayRead(vecs[i],(const PetscScalar **)&petscvecarray);
177:         VecHYPRE_ParVectorReplacePointer(jac->hmnull[i],petscvecarray,jac->hmnull_hypre_data_array[i]);
178:         VecRestoreArrayRead(vecs[i],(const PetscScalar **)&petscvecarray);
179:         PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->hmnull[i],(void**)&jac->phmnull[i]));
180:       }
181:       if (has_const) {
182:         MatCreateVecs(pc->pmat,&jac->hmnull_constant,NULL);
183:         VecSet(jac->hmnull_constant,1);
184:         VecNormalize(jac->hmnull_constant,NULL);
185:         VecHYPRE_IJVectorCreate(jac->hmnull_constant,&jac->hmnull[nvec]);
186:         VecGetArrayRead(jac->hmnull_constant,(const PetscScalar **)&petscvecarray);
187:         VecHYPRE_ParVectorReplacePointer(jac->hmnull[nvec],petscvecarray,jac->hmnull_hypre_data_array[nvec]);
188:         VecRestoreArrayRead(jac->hmnull_constant,(const PetscScalar **)&petscvecarray);
189:         PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->hmnull[nvec],(void**)&jac->phmnull[nvec]));
190:         nvec++;
191:       }
192:       PetscStackCallStandard(HYPRE_BoomerAMGSetInterpVectors,(jac->hsolver,nvec,jac->phmnull));
193:       jac->n_hmnull = nvec;
194:     }
195:   }

197:   /* special case for AMS */
198:   if (jac->setup == HYPRE_AMSSetup) {
199:     Mat_HYPRE          *hm;
200:     HYPRE_ParCSRMatrix parcsr;
201:     if (!jac->coords[0] && !jac->constants[0] && !(jac->ND_PiFull || (jac->ND_Pi[0] && jac->ND_Pi[1]))) {
202:       SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE AMS preconditioner needs either the coordinate vectors via PCSetCoordinates() or the edge constant vectors via PCHYPRESetEdgeConstantVectors() or the interpolation matrix via PCHYPRESetInterpolations");
203:     }
204:     if (jac->dim) {
205:       PetscStackCallStandard(HYPRE_AMSSetDimension,(jac->hsolver,jac->dim));
206:     }
207:     if (jac->constants[0]) {
208:       HYPRE_ParVector ozz,zoz,zzo = NULL;
209:       PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->constants[0],(void**)(&ozz)));
210:       PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->constants[1],(void**)(&zoz)));
211:       if (jac->constants[2]) {
212:         PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->constants[2],(void**)(&zzo)));
213:       }
214:       PetscStackCallStandard(HYPRE_AMSSetEdgeConstantVectors,(jac->hsolver,ozz,zoz,zzo));
215:     }
216:     if (jac->coords[0]) {
217:       HYPRE_ParVector coords[3];
218:       coords[0] = NULL;
219:       coords[1] = NULL;
220:       coords[2] = NULL;
221:       if (jac->coords[0]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[0],(void**)(&coords[0])));
222:       if (jac->coords[1]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[1],(void**)(&coords[1])));
223:       if (jac->coords[2]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[2],(void**)(&coords[2])));
224:       PetscStackCallStandard(HYPRE_AMSSetCoordinateVectors,(jac->hsolver,coords[0],coords[1],coords[2]));
225:     }
226:     if (!jac->G) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE AMS preconditioner needs the discrete gradient operator via PCHYPRESetDiscreteGradient");
227:     hm = (Mat_HYPRE*)(jac->G->data);
228:     PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&parcsr)));
229:     PetscStackCallStandard(HYPRE_AMSSetDiscreteGradient,(jac->hsolver,parcsr));
230:     if (jac->alpha_Poisson) {
231:       hm = (Mat_HYPRE*)(jac->alpha_Poisson->data);
232:       PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&parcsr)));
233:       PetscStackCallStandard(HYPRE_AMSSetAlphaPoissonMatrix,(jac->hsolver,parcsr));
234:     }
235:     if (jac->ams_beta_is_zero) {
236:       PetscStackCallStandard(HYPRE_AMSSetBetaPoissonMatrix,(jac->hsolver,NULL));
237:     } else if (jac->beta_Poisson) {
238:       hm = (Mat_HYPRE*)(jac->beta_Poisson->data);
239:       PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&parcsr)));
240:       PetscStackCallStandard(HYPRE_AMSSetBetaPoissonMatrix,(jac->hsolver,parcsr));
241:     }
242:     if (jac->ND_PiFull || (jac->ND_Pi[0] && jac->ND_Pi[1])) {
243:       PetscInt           i;
244:       HYPRE_ParCSRMatrix nd_parcsrfull, nd_parcsr[3];
245:       if (jac->ND_PiFull) {
246:         hm = (Mat_HYPRE*)(jac->ND_PiFull->data);
247:         PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&nd_parcsrfull)));
248:       } else {
249:         nd_parcsrfull = NULL;
250:       }
251:       for (i=0;i<3;++i) {
252:         if (jac->ND_Pi[i]) {
253:           hm = (Mat_HYPRE*)(jac->ND_Pi[i]->data);
254:           PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&nd_parcsr[i])));
255:         } else {
256:           nd_parcsr[i] = NULL;
257:         }
258:       }
259:       PetscStackCallStandard(HYPRE_AMSSetInterpolations,(jac->hsolver,nd_parcsrfull,nd_parcsr[0],nd_parcsr[1],nd_parcsr[2]));
260:     }
261:   }
262:   /* special case for ADS */
263:   if (jac->setup == HYPRE_ADSSetup) {
264:     Mat_HYPRE          *hm;
265:     HYPRE_ParCSRMatrix parcsr;
266:     if (!jac->coords[0] && !((jac->RT_PiFull || (jac->RT_Pi[0] && jac->RT_Pi[1])) && (jac->ND_PiFull || (jac->ND_Pi[0] && jac->ND_Pi[1])))) {
267:       SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE ADS preconditioner needs either the coordinate vectors via PCSetCoordinates() or the interpolation matrices via PCHYPRESetInterpolations");
268:     }
269:     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");
270:     if (!jac->G) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE ADS preconditioner needs the discrete gradient operator via PCHYPRESetDiscreteGradient");
271:     if (!jac->C) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"HYPRE ADS preconditioner needs the discrete curl operator via PCHYPRESetDiscreteGradient");
272:     if (jac->coords[0]) {
273:       HYPRE_ParVector coords[3];
274:       coords[0] = NULL;
275:       coords[1] = NULL;
276:       coords[2] = NULL;
277:       if (jac->coords[0]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[0],(void**)(&coords[0])));
278:       if (jac->coords[1]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[1],(void**)(&coords[1])));
279:       if (jac->coords[2]) PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->coords[2],(void**)(&coords[2])));
280:       PetscStackCallStandard(HYPRE_ADSSetCoordinateVectors,(jac->hsolver,coords[0],coords[1],coords[2]));
281:     }
282:     hm = (Mat_HYPRE*)(jac->G->data);
283:     PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&parcsr)));
284:     PetscStackCallStandard(HYPRE_ADSSetDiscreteGradient,(jac->hsolver,parcsr));
285:     hm = (Mat_HYPRE*)(jac->C->data);
286:     PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&parcsr)));
287:     PetscStackCallStandard(HYPRE_ADSSetDiscreteCurl,(jac->hsolver,parcsr));
288:     if ((jac->RT_PiFull || (jac->RT_Pi[0] && jac->RT_Pi[1])) && (jac->ND_PiFull || (jac->ND_Pi[0] && jac->ND_Pi[1]))) {
289:       PetscInt           i;
290:       HYPRE_ParCSRMatrix rt_parcsrfull, rt_parcsr[3];
291:       HYPRE_ParCSRMatrix nd_parcsrfull, nd_parcsr[3];
292:       if (jac->RT_PiFull) {
293:         hm = (Mat_HYPRE*)(jac->RT_PiFull->data);
294:         PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&rt_parcsrfull)));
295:       } else {
296:         rt_parcsrfull = NULL;
297:       }
298:       for (i=0;i<3;++i) {
299:         if (jac->RT_Pi[i]) {
300:           hm = (Mat_HYPRE*)(jac->RT_Pi[i]->data);
301:           PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&rt_parcsr[i])));
302:         } else {
303:           rt_parcsr[i] = NULL;
304:         }
305:       }
306:       if (jac->ND_PiFull) {
307:         hm = (Mat_HYPRE*)(jac->ND_PiFull->data);
308:         PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&nd_parcsrfull)));
309:       } else {
310:         nd_parcsrfull = NULL;
311:       }
312:       for (i=0;i<3;++i) {
313:         if (jac->ND_Pi[i]) {
314:           hm = (Mat_HYPRE*)(jac->ND_Pi[i]->data);
315:           PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hm->ij,(void**)(&nd_parcsr[i])));
316:         } else {
317:           nd_parcsr[i] = NULL;
318:         }
319:       }
320:       PetscStackCallStandard(HYPRE_ADSSetInterpolations,(jac->hsolver,rt_parcsrfull,rt_parcsr[0],rt_parcsr[1],rt_parcsr[2],nd_parcsrfull,nd_parcsr[0],nd_parcsr[1],nd_parcsr[2]));
321:     }
322:   }
323:   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hjac->ij,(void**)&hmat));
324:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hjac->b,(void**)&bv));
325:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hjac->x,(void**)&xv));
326:   PetscStackCallStandard(jac->setup,(jac->hsolver,hmat,bv,xv));
327:   return(0);
328: }

330: static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
331: {
332:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
333:   Mat_HYPRE          *hjac = (Mat_HYPRE*)(jac->hpmat->data);
334:   PetscErrorCode     ierr;
335:   HYPRE_ParCSRMatrix hmat;
336:   PetscScalar        *xv;
337:   const PetscScalar  *bv,*sbv;
338:   HYPRE_ParVector    jbv,jxv;
339:   PetscScalar        *sxv;
340:   PetscInt           hierr;

343:   PetscCitationsRegister(hypreCitation,&cite);
344:   if (!jac->applyrichardson) {VecSet(x,0.0);}
345:   VecGetArrayRead(b,&bv);
346:   VecGetArray(x,&xv);
347:   VecHYPRE_ParVectorReplacePointer(hjac->b,(PetscScalar*)bv,sbv);
348:   VecHYPRE_ParVectorReplacePointer(hjac->x,xv,sxv);

350:   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hjac->ij,(void**)&hmat));
351:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hjac->b,(void**)&jbv));
352:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hjac->x,(void**)&jxv));
353:   PetscStackCall("Hypre solve",h(*jac->solve)(jac->hsolver,hmat,jbv,jxv);
354:   if (hierr && hierr != HYPRE_ERROR_CONV) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
355:   if (hierr) hypre__global_error = 0;);

357:   if (jac->setup == HYPRE_AMSSetup && jac->ams_beta_is_zero_part) {
358:     PetscStackCallStandard(HYPRE_AMSProjectOutGradients,(jac->hsolver,jxv));
359:   }
360:   VecHYPRE_ParVectorReplacePointer(hjac->b,(PetscScalar*)sbv,bv);
361:   VecHYPRE_ParVectorReplacePointer(hjac->x,sxv,xv);
362:   VecRestoreArray(x,&xv);
363:   VecRestoreArrayRead(b,&bv);
364:   return(0);
365: }

367: static PetscErrorCode PCReset_HYPRE(PC pc)
368: {
369:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;

373:   MatDestroy(&jac->hpmat);
374:   MatDestroy(&jac->G);
375:   MatDestroy(&jac->C);
376:   MatDestroy(&jac->alpha_Poisson);
377:   MatDestroy(&jac->beta_Poisson);
378:   MatDestroy(&jac->RT_PiFull);
379:   MatDestroy(&jac->RT_Pi[0]);
380:   MatDestroy(&jac->RT_Pi[1]);
381:   MatDestroy(&jac->RT_Pi[2]);
382:   MatDestroy(&jac->ND_PiFull);
383:   MatDestroy(&jac->ND_Pi[0]);
384:   MatDestroy(&jac->ND_Pi[1]);
385:   MatDestroy(&jac->ND_Pi[2]);
386:   if (jac->coords[0]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[0])); jac->coords[0] = NULL;
387:   if (jac->coords[1]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[1])); jac->coords[1] = NULL;
388:   if (jac->coords[2]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[2])); jac->coords[2] = NULL;
389:   if (jac->constants[0]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[0])); jac->constants[0] = NULL;
390:   if (jac->constants[1]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[1])); jac->constants[1] = NULL;
391:   if (jac->constants[2]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[2])); jac->constants[2] = NULL;
392:   if (jac->n_hmnull && jac->hmnull) {
393:     PetscInt                 i;
394:     PETSC_UNUSED PetscScalar *petscvecarray;

396:     for (i=0; i<jac->n_hmnull; i++) {
397:       VecHYPRE_ParVectorReplacePointer(jac->hmnull[i],jac->hmnull_hypre_data_array[i],petscvecarray);
398:       PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->hmnull[i]));
399:     }
400:     PetscFree(jac->hmnull);
401:     PetscFree(jac->hmnull_hypre_data_array);
402:     PetscFree(jac->phmnull);
403:     VecDestroy(&jac->hmnull_constant);
404:   }
405:   jac->ams_beta_is_zero = PETSC_FALSE;
406:   jac->dim = 0;
407:   return(0);
408: }

410: static PetscErrorCode PCDestroy_HYPRE(PC pc)
411: {
412:   PC_HYPRE                 *jac = (PC_HYPRE*)pc->data;
413:   PetscErrorCode           ierr;

416:   PCReset_HYPRE(pc);
417:   if (jac->destroy) PetscStackCallStandard(jac->destroy,(jac->hsolver));
418:   PetscFree(jac->hypre_type);
419:   if (jac->comm_hypre != MPI_COMM_NULL) { MPI_Comm_free(&(jac->comm_hypre));}
420:   PetscFree(pc->data);

422:   PetscObjectChangeTypeName((PetscObject)pc,0);
423:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetType_C",NULL);
424:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPREGetType_C",NULL);
425:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetCoordinates_C",NULL);
426:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteGradient_C",NULL);
427:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteCurl_C",NULL);
428:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetInterpolations_C",NULL);
429:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetConstantEdgeVectors_C",NULL);
430:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetPoissonMatrix_C",NULL);
431:   return(0);
432: }

434: /* --------------------------------------------------------------------------------------------*/
435: static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PetscOptionItems *PetscOptionsObject,PC pc)
436: {
437:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
439:   PetscBool      flag;

442:   PetscOptionsHead(PetscOptionsObject,"HYPRE Pilut Options");
443:   PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);
444:   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetMaxIter,(jac->hsolver,jac->maxiter));
445:   PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);
446:   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetDropTolerance,(jac->hsolver,jac->tol));
447:   PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);
448:   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetFactorRowSize,(jac->hsolver,jac->factorrowsize));
449:   PetscOptionsTail();
450:   return(0);
451: }

453: static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
454: {
455:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
457:   PetscBool      iascii;

460:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
461:   if (iascii) {
462:     PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");
463:     if (jac->maxiter != PETSC_DEFAULT) {
464:       PetscViewerASCIIPrintf(viewer,"    maximum number of iterations %d\n",jac->maxiter);
465:     } else {
466:       PetscViewerASCIIPrintf(viewer,"    default maximum number of iterations \n");
467:     }
468:     if (jac->tol != PETSC_DEFAULT) {
469:       PetscViewerASCIIPrintf(viewer,"    drop tolerance %g\n",(double)jac->tol);
470:     } else {
471:       PetscViewerASCIIPrintf(viewer,"    default drop tolerance \n");
472:     }
473:     if (jac->factorrowsize != PETSC_DEFAULT) {
474:       PetscViewerASCIIPrintf(viewer,"    factor row size %d\n",jac->factorrowsize);
475:     } else {
476:       PetscViewerASCIIPrintf(viewer,"    default factor row size \n");
477:     }
478:   }
479:   return(0);
480: }

482: /* --------------------------------------------------------------------------------------------*/
483: static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PetscOptionItems *PetscOptionsObject,PC pc)
484: {
485:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
487:   PetscBool      flag;

490:   PetscOptionsHead(PetscOptionsObject,"HYPRE Euclid Options");
491:   PetscOptionsInt("-pc_hypre_euclid_level","Factorization levels","None",jac->eu_level,&jac->eu_level,&flag);
492:   if (flag) PetscStackCallStandard(HYPRE_EuclidSetLevel,(jac->hsolver,jac->eu_level));
493:   PetscOptionsTail();
494:   return(0);
495: }

497: static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
498: {
499:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
501:   PetscBool      iascii;

504:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
505:   if (iascii) {
506:     PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n");
507:     if (jac->eu_level != PETSC_DEFAULT) {
508:       PetscViewerASCIIPrintf(viewer,"    factorization levels %d\n",jac->eu_level);
509:     } else {
510:       PetscViewerASCIIPrintf(viewer,"    default factorization levels \n");
511:     }
512:   }
513:   return(0);
514: }

516: /* --------------------------------------------------------------------------------------------*/

518: static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
519: {
520:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
521:   Mat_HYPRE          *hjac = (Mat_HYPRE*)(jac->hpmat->data);
522:   PetscErrorCode     ierr;
523:   HYPRE_ParCSRMatrix hmat;
524:   PetscScalar        *xv;
525:   const PetscScalar  *bv;
526:   HYPRE_ParVector    jbv,jxv;
527:   PetscScalar        *sbv,*sxv;
528:   PetscInt           hierr;

531:   PetscCitationsRegister(hypreCitation,&cite);
532:   VecSet(x,0.0);
533:   VecGetArrayRead(b,&bv);
534:   VecGetArray(x,&xv);
535:   VecHYPRE_ParVectorReplacePointer(hjac->b,(PetscScalar*)bv,sbv);
536:   VecHYPRE_ParVectorReplacePointer(hjac->x,xv,sxv);

538:   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hjac->ij,(void**)&hmat));
539:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hjac->b,(void**)&jbv));
540:   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hjac->x,(void**)&jxv));

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

547:   VecHYPRE_ParVectorReplacePointer(hjac->b,sbv,bv);
548:   VecHYPRE_ParVectorReplacePointer(hjac->x,sxv,xv);
549:   VecRestoreArray(x,&xv);
550:   VecRestoreArrayRead(b,&bv);
551:   return(0);
552: }

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

557: static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
558: static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"};
559: static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
560: /* The following corresponds to HYPRE_BoomerAMGSetRelaxType which has many missing numbers in the enum */
561: static const char *HYPREBoomerAMGSmoothType[]   = {"Schwarz-smoothers","Pilut","ParaSails","Euclid"};
562: static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","seqboundary-Gauss-Seidel","SOR/Jacobi","backward-SOR/Jacobi",
563:                                                   "" /* [5] hybrid chaotic Gauss-Seidel (works only with OpenMP) */,"symmetric-SOR/Jacobi",
564:                                                   "" /* 7 */,"l1scaled-SOR/Jacobi","Gaussian-elimination",
565:                                                   "" /* 10 */, "" /* 11 */, "" /* 12 */, "l1-Gauss-Seidel" /* nonsymmetric */, "backward-l1-Gauss-Seidel" /* nonsymmetric */,
566:                                                   "CG" /* non-stationary */,"Chebyshev","FCF-Jacobi","l1scaled-Jacobi"};
567: static const char *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i",
568:                                                   "ext+i-cc", "standard", "standard-wts", "block", "block-wtd", "FF", "FF1"};
569: static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PetscOptionItems *PetscOptionsObject,PC pc)
570: {
571:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
573:   PetscInt       bs,n,indx,level;
574:   PetscBool      flg, tmp_truth;
575:   double         tmpdbl, twodbl[2];

578:   PetscOptionsHead(PetscOptionsObject,"HYPRE BoomerAMG Options");
579:   PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType+1,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);
580:   if (flg) {
581:     jac->cycletype = indx+1;
582:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
583:   }
584:   PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);
585:   if (flg) {
586:     if (jac->maxlevels < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
587:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
588:   }
589:   PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);
590:   if (flg) {
591:     if (jac->maxiter < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
592:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
593:   }
594:   PetscOptionsScalar("-pc_hypre_boomeramg_tol","Convergence tolerance PER hypre call (0.0 = use a fixed number of iterations)","None",jac->tol,&jac->tol,&flg);
595:   if (flg) {
596:     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);
597:     PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
598:   }
599:   bs = 1;
600:   if (pc->pmat) {
601:     MatGetBlockSize(pc->pmat,&bs);
602:   }
603:   PetscOptionsInt("-pc_hypre_boomeramg_numfunctions","Number of functions","HYPRE_BoomerAMGSetNumFunctions",bs,&bs,&flg);
604:   if (flg) {
605:     PetscStackCallStandard(HYPRE_BoomerAMGSetNumFunctions,(jac->hsolver,bs));
606:   }

608:   PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);
609:   if (flg) {
610:     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);
611:     PetscStackCallStandard(HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
612:   }

614:   PetscOptionsInt("-pc_hypre_boomeramg_P_max","Max elements per row for interpolation operator (0=unlimited)","None",jac->pmax,&jac->pmax,&flg);
615:   if (flg) {
616:     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);
617:     PetscStackCallStandard(HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
618:   }

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

624:     PetscStackCallStandard(HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
625:   }

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

631:     PetscStackCallStandard(HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
632:   }


635:   PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);
636:   if (flg) {
637:     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);
638:     PetscStackCallStandard(HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
639:   }
640:   PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);
641:   if (flg) {
642:     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);
643:     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);
644:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
645:   }

647:   /* Grid sweeps */
648:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_all","Number of sweeps for the up and down grid levels","None",jac->gridsweeps[0],&indx,&flg);
649:   if (flg) {
650:     PetscStackCallStandard(HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver,indx));
651:     /* modify the jac structure so we can view the updated options with PC_View */
652:     jac->gridsweeps[0] = indx;
653:     jac->gridsweeps[1] = indx;
654:     /*defaults coarse to 1 */
655:     jac->gridsweeps[2] = 1;
656:   }
657:   PetscOptionsInt("-pc_hypre_boomeramg_nodal_coarsen","Use a nodal based coarsening 1-6","HYPRE_BoomerAMGSetNodal",jac->nodal_coarsening,&jac->nodal_coarsening,&flg);
658:   if (flg) {
659:     PetscStackCallStandard(HYPRE_BoomerAMGSetNodal,(jac->hsolver,jac->nodal_coarsening));
660:   }
661:   PetscOptionsInt("-pc_hypre_boomeramg_nodal_coarsen_diag","Diagonal in strength matrix for nodal based coarsening 0-2","HYPRE_BoomerAMGSetNodalDiag",jac->nodal_coarsening_diag,&jac->nodal_coarsening_diag,&flg);
662:   if (flg) {
663:     PetscStackCallStandard(HYPRE_BoomerAMGSetNodalDiag,(jac->hsolver,jac->nodal_coarsening_diag));
664:   }
665:   PetscOptionsInt("-pc_hypre_boomeramg_vec_interp_variant","Variant of algorithm 1-3","HYPRE_BoomerAMGSetInterpVecVariant",jac->vec_interp_variant, &jac->vec_interp_variant,&flg);
666:   if (flg) {
667:     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpVecVariant,(jac->hsolver,jac->vec_interp_variant));
668:   }
669:   PetscOptionsInt("-pc_hypre_boomeramg_vec_interp_qmax","Max elements per row for each Q","HYPRE_BoomerAMGSetInterpVecQMax",jac->vec_interp_qmax, &jac->vec_interp_qmax,&flg);
670:   if (flg) {
671:     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpVecQMax,(jac->hsolver,jac->vec_interp_qmax));
672:   }
673:   PetscOptionsBool("-pc_hypre_boomeramg_vec_interp_smooth","Whether to smooth the interpolation vectors","HYPRE_BoomerAMGSetSmoothInterpVectors",jac->vec_interp_smooth, &jac->vec_interp_smooth,&flg);
674:   if (flg) {
675:     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothInterpVectors,(jac->hsolver,jac->vec_interp_smooth));
676:   }
677:   PetscOptionsInt("-pc_hypre_boomeramg_interp_refine","Preprocess the interpolation matrix through iterative weight refinement","HYPRE_BoomerAMGSetInterpRefine",jac->interp_refine, &jac->interp_refine,&flg);
678:   if (flg) {
679:     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpRefine,(jac->hsolver,jac->interp_refine));
680:   }
681:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx,&flg);
682:   if (flg) {
683:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 1));
684:     jac->gridsweeps[0] = indx;
685:   }
686:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);
687:   if (flg) {
688:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 2));
689:     jac->gridsweeps[1] = indx;
690:   }
691:   PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);
692:   if (flg) {
693:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 3));
694:     jac->gridsweeps[2] = indx;
695:   }

697:   /* Smooth type */
698:   PetscOptionsEList("-pc_hypre_boomeramg_smooth_type","Enable more complex smoothers","None",HYPREBoomerAMGSmoothType,ALEN(HYPREBoomerAMGSmoothType),HYPREBoomerAMGSmoothType[0],&indx,&flg);
699:   if (flg) {
700:     jac->smoothtype = indx;
701:     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothType,(jac->hsolver,indx+6));
702:     jac->smoothnumlevels = 25;
703:     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothNumLevels,(jac->hsolver,25));
704:   }

706:   /* Number of smoothing levels */
707:   PetscOptionsInt("-pc_hypre_boomeramg_smooth_num_levels","Number of levels on which more complex smoothers are used","None",25,&indx,&flg);
708:   if (flg && (jac->smoothtype != -1)) {
709:     jac->smoothnumlevels = indx;
710:     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothNumLevels,(jac->hsolver,indx));
711:   }

713:   /* Number of levels for ILU(k) for Euclid */
714:   PetscOptionsInt("-pc_hypre_boomeramg_eu_level","Number of levels for ILU(k) in Euclid smoother","None",0,&indx,&flg);
715:   if (flg && (jac->smoothtype == 3)) {
716:     jac->eu_level = indx;
717:     PetscStackCallStandard(HYPRE_BoomerAMGSetEuLevel,(jac->hsolver,indx));
718:   }

720:   /* Filter for ILU(k) for Euclid */
721:   double droptolerance;
722:   PetscOptionsScalar("-pc_hypre_boomeramg_eu_droptolerance","Drop tolerance for ILU(k) in Euclid smoother","None",0,&droptolerance,&flg);
723:   if (flg && (jac->smoothtype == 3)) {
724:     jac->eu_droptolerance = droptolerance;
725:     PetscStackCallStandard(HYPRE_BoomerAMGSetEuLevel,(jac->hsolver,droptolerance));
726:   }

728:   /* Use Block Jacobi ILUT for Euclid */
729:   PetscOptionsBool("-pc_hypre_boomeramg_eu_bj", "Use Block Jacobi for ILU in Euclid smoother?", "None", PETSC_FALSE, &tmp_truth, &flg);
730:   if (flg && (jac->smoothtype == 3)) {
731:     jac->eu_bj = tmp_truth;
732:     PetscStackCallStandard(HYPRE_BoomerAMGSetEuBJ,(jac->hsolver,jac->eu_bj));
733:   }

735:   /* Relax type */
736:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_all","Relax type for the up and down cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);
737:   if (flg) {
738:     jac->relaxtype[0] = jac->relaxtype[1]  = indx;
739:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, indx));
740:     /* by default, coarse type set to 9 */
741:     jac->relaxtype[2] = 9;
742:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, 9, 3));
743:   }
744:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);
745:   if (flg) {
746:     jac->relaxtype[0] = indx;
747:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 1));
748:   }
749:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);
750:   if (flg) {
751:     jac->relaxtype[1] = indx;
752:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 2));
753:   }
754:   PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[9],&indx,&flg);
755:   if (flg) {
756:     jac->relaxtype[2] = indx;
757:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 3));
758:   }

760:   /* Relaxation Weight */
761:   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);
762:   if (flg) {
763:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxWt,(jac->hsolver,tmpdbl));
764:     jac->relaxweight = tmpdbl;
765:   }

767:   n         = 2;
768:   twodbl[0] = twodbl[1] = 1.0;
769:   PetscOptionsRealArray("-pc_hypre_boomeramg_relax_weight_level","Set the relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);
770:   if (flg) {
771:     if (n == 2) {
772:       indx =  (int)PetscAbsReal(twodbl[1]);
773:       PetscStackCallStandard(HYPRE_BoomerAMGSetLevelRelaxWt,(jac->hsolver,twodbl[0],indx));
774:     } 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);
775:   }

777:   /* Outer relaxation Weight */
778:   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);
779:   if (flg) {
780:     PetscStackCallStandard(HYPRE_BoomerAMGSetOuterWt,(jac->hsolver, tmpdbl));
781:     jac->outerrelaxweight = tmpdbl;
782:   }

784:   n         = 2;
785:   twodbl[0] = twodbl[1] = 1.0;
786:   PetscOptionsRealArray("-pc_hypre_boomeramg_outer_relax_weight_level","Set the outer relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);
787:   if (flg) {
788:     if (n == 2) {
789:       indx =  (int)PetscAbsReal(twodbl[1]);
790:       PetscStackCallStandard(HYPRE_BoomerAMGSetLevelOuterWt,(jac->hsolver, twodbl[0], indx));
791:     } 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);
792:   }

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

797:   if (flg && tmp_truth) {
798:     jac->relaxorder = 0;
799:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
800:   }
801:   PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,ALEN(HYPREBoomerAMGMeasureType),HYPREBoomerAMGMeasureType[0],&indx,&flg);
802:   if (flg) {
803:     jac->measuretype = indx;
804:     PetscStackCallStandard(HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
805:   }
806:   /* update list length 3/07 */
807:   PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,ALEN(HYPREBoomerAMGCoarsenType),HYPREBoomerAMGCoarsenType[6],&indx,&flg);
808:   if (flg) {
809:     jac->coarsentype = indx;
810:     PetscStackCallStandard(HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
811:   }

813:   /* new 3/07 */
814:   PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,ALEN(HYPREBoomerAMGInterpType),HYPREBoomerAMGInterpType[0],&indx,&flg);
815:   if (flg) {
816:     jac->interptype = indx;
817:     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
818:   }

820:   PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);
821:   if (flg) {
822:     level = 3;
823:     PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,NULL);

825:     jac->printstatistics = PETSC_TRUE;
826:     PetscStackCallStandard(HYPRE_BoomerAMGSetPrintLevel,(jac->hsolver,level));
827:   }

829:   PetscOptionsName("-pc_hypre_boomeramg_print_debug","Print debug information","None",&flg);
830:   if (flg) {
831:     level = 3;
832:     PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,NULL);

834:     jac->printstatistics = PETSC_TRUE;
835:     PetscStackCallStandard(HYPRE_BoomerAMGSetDebugFlag,(jac->hsolver,level));
836:   }

838:   PetscOptionsBool("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);
839:   if (flg && tmp_truth) {
840:     PetscInt tmp_int;
841:     PetscOptionsInt("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);
842:     if (flg) jac->nodal_relax_levels = tmp_int;
843:     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothType,(jac->hsolver,6));
844:     PetscStackCallStandard(HYPRE_BoomerAMGSetDomainType,(jac->hsolver,1));
845:     PetscStackCallStandard(HYPRE_BoomerAMGSetOverlap,(jac->hsolver,0));
846:     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothNumLevels,(jac->hsolver,jac->nodal_relax_levels));
847:   }

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

853: 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)
854: {
855:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
857:   PetscInt       oits;

860:   PetscCitationsRegister(hypreCitation,&cite);
861:   PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,its*jac->maxiter));
862:   PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,rtol));
863:   jac->applyrichardson = PETSC_TRUE;
864:   PCApply_HYPRE(pc,b,y);
865:   jac->applyrichardson = PETSC_FALSE;
866:   PetscStackCallStandard(HYPRE_BoomerAMGGetNumIterations,(jac->hsolver,(HYPRE_Int *)&oits));
867:   *outits = oits;
868:   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
869:   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
870:   PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
871:   PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
872:   return(0);
873: }


876: static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
877: {
878:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
880:   PetscBool      iascii;

883:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
884:   if (iascii) {
885:     PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");
886:     PetscViewerASCIIPrintf(viewer,"    Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);
887:     PetscViewerASCIIPrintf(viewer,"    Maximum number of levels %D\n",jac->maxlevels);
888:     PetscViewerASCIIPrintf(viewer,"    Maximum number of iterations PER hypre call %D\n",jac->maxiter);
889:     PetscViewerASCIIPrintf(viewer,"    Convergence tolerance PER hypre call %g\n",(double)jac->tol);
890:     PetscViewerASCIIPrintf(viewer,"    Threshold for strong coupling %g\n",(double)jac->strongthreshold);
891:     PetscViewerASCIIPrintf(viewer,"    Interpolation truncation factor %g\n",(double)jac->truncfactor);
892:     PetscViewerASCIIPrintf(viewer,"    Interpolation: max elements per row %D\n",jac->pmax);
893:     if (jac->interp_refine) {
894:       PetscViewerASCIIPrintf(viewer,"    Interpolation: number of steps of weighted refinement %D\n",jac->interp_refine);
895:     }
896:     PetscViewerASCIIPrintf(viewer,"    Number of levels of aggressive coarsening %D\n",jac->agg_nl);
897:     PetscViewerASCIIPrintf(viewer,"    Number of paths for aggressive coarsening %D\n",jac->agg_num_paths);

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

901:     PetscViewerASCIIPrintf(viewer,"    Sweeps down         %D\n",jac->gridsweeps[0]);
902:     PetscViewerASCIIPrintf(viewer,"    Sweeps up           %D\n",jac->gridsweeps[1]);
903:     PetscViewerASCIIPrintf(viewer,"    Sweeps on coarse    %D\n",jac->gridsweeps[2]);

905:     PetscViewerASCIIPrintf(viewer,"    Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);
906:     PetscViewerASCIIPrintf(viewer,"    Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);
907:     PetscViewerASCIIPrintf(viewer,"    Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);

909:     PetscViewerASCIIPrintf(viewer,"    Relax weight  (all)      %g\n",(double)jac->relaxweight);
910:     PetscViewerASCIIPrintf(viewer,"    Outer relax weight (all) %g\n",(double)jac->outerrelaxweight);

912:     if (jac->relaxorder) {
913:       PetscViewerASCIIPrintf(viewer,"    Using CF-relaxation\n");
914:     } else {
915:       PetscViewerASCIIPrintf(viewer,"    Not using CF-relaxation\n");
916:     }
917:     if (jac->smoothtype!=-1) {
918:       PetscViewerASCIIPrintf(viewer,"    Smooth type          %s\n",HYPREBoomerAMGSmoothType[jac->smoothtype]);
919:       PetscViewerASCIIPrintf(viewer,"    Smooth num levels    %D\n",jac->smoothnumlevels);
920:     } else {
921:       PetscViewerASCIIPrintf(viewer,"    Not using more complex smoothers.\n");
922:     }
923:     if (jac->smoothtype==3) {
924:       PetscViewerASCIIPrintf(viewer,"    Euclid ILU(k) levels %D\n",jac->eu_level);
925:       PetscViewerASCIIPrintf(viewer,"    Euclid ILU(k) drop tolerance %g\n",(double)jac->eu_droptolerance);
926:       PetscViewerASCIIPrintf(viewer,"    Euclid ILU use Block-Jacobi? %D\n",jac->eu_bj);
927:     }
928:     PetscViewerASCIIPrintf(viewer,"    Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);
929:     PetscViewerASCIIPrintf(viewer,"    Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);
930:     PetscViewerASCIIPrintf(viewer,"    Interpolation type  %s\n",HYPREBoomerAMGInterpType[jac->interptype]);
931:     if (jac->nodal_coarsening) {
932:       PetscViewerASCIIPrintf(viewer,"    Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal() %D\n",jac->nodal_coarsening);
933:     }
934:     if (jac->vec_interp_variant) {
935:       PetscViewerASCIIPrintf(viewer,"    HYPRE_BoomerAMGSetInterpVecVariant() %D\n",jac->vec_interp_variant);
936:       PetscViewerASCIIPrintf(viewer,"    HYPRE_BoomerAMGSetInterpVecQMax() %D\n",jac->vec_interp_qmax);
937:       PetscViewerASCIIPrintf(viewer,"    HYPRE_BoomerAMGSetSmoothInterpVectors() %d\n",jac->vec_interp_smooth);
938:     }
939:     if (jac->nodal_relax) {
940:       PetscViewerASCIIPrintf(viewer,"    Using nodal relaxation via Schwarz smoothing on levels %D\n",jac->nodal_relax_levels);
941:     }
942:   }
943:   return(0);
944: }

946: /* --------------------------------------------------------------------------------------------*/
947: static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PetscOptionItems *PetscOptionsObject,PC pc)
948: {
949:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
951:   PetscInt       indx;
952:   PetscBool      flag;
953:   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};

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

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

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

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

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

973:   PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,ALEN(symtlist),symtlist[0],&indx,&flag);
974:   if (flag) {
975:     jac->symt = indx;
976:     PetscStackCallStandard(HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
977:   }

979:   PetscOptionsTail();
980:   return(0);
981: }

983: static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
984: {
985:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
987:   PetscBool      iascii;
988:   const char     *symt = 0;;

991:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
992:   if (iascii) {
993:     PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");
994:     PetscViewerASCIIPrintf(viewer,"    nlevels %d\n",jac->nlevels);
995:     PetscViewerASCIIPrintf(viewer,"    threshold %g\n",(double)jac->threshhold);
996:     PetscViewerASCIIPrintf(viewer,"    filter %g\n",(double)jac->filter);
997:     PetscViewerASCIIPrintf(viewer,"    load balance %g\n",(double)jac->loadbal);
998:     PetscViewerASCIIPrintf(viewer,"    reuse nonzero structure %s\n",PetscBools[jac->ruse]);
999:     PetscViewerASCIIPrintf(viewer,"    print info to screen %s\n",PetscBools[jac->logging]);
1000:     if (!jac->symt) symt = "nonsymmetric matrix and preconditioner";
1001:     else if (jac->symt == 1) symt = "SPD matrix and preconditioner";
1002:     else if (jac->symt == 2) symt = "nonsymmetric matrix but SPD preconditioner";
1003:     else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
1004:     PetscViewerASCIIPrintf(viewer,"    %s\n",symt);
1005:   }
1006:   return(0);
1007: }
1008: /* --------------------------------------------------------------------------------------------*/
1009: static PetscErrorCode PCSetFromOptions_HYPRE_AMS(PetscOptionItems *PetscOptionsObject,PC pc)
1010: {
1011:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1013:   PetscInt       n;
1014:   PetscBool      flag,flag2,flag3,flag4;

1017:   PetscOptionsHead(PetscOptionsObject,"HYPRE AMS Options");
1018:   PetscOptionsInt("-pc_hypre_ams_print_level","Debugging output level for AMS","None",jac->as_print,&jac->as_print,&flag);
1019:   if (flag) PetscStackCallStandard(HYPRE_AMSSetPrintLevel,(jac->hsolver,jac->as_print));
1020:   PetscOptionsInt("-pc_hypre_ams_max_iter","Maximum number of AMS multigrid iterations within PCApply","None",jac->as_max_iter,&jac->as_max_iter,&flag);
1021:   if (flag) PetscStackCallStandard(HYPRE_AMSSetMaxIter,(jac->hsolver,jac->as_max_iter));
1022:   PetscOptionsInt("-pc_hypre_ams_cycle_type","Cycle type for AMS multigrid","None",jac->ams_cycle_type,&jac->ams_cycle_type,&flag);
1023:   if (flag) PetscStackCallStandard(HYPRE_AMSSetCycleType,(jac->hsolver,jac->ams_cycle_type));
1024:   PetscOptionsReal("-pc_hypre_ams_tol","Error tolerance for AMS multigrid","None",jac->as_tol,&jac->as_tol,&flag);
1025:   if (flag) PetscStackCallStandard(HYPRE_AMSSetTol,(jac->hsolver,jac->as_tol));
1026:   PetscOptionsInt("-pc_hypre_ams_relax_type","Relaxation type for AMS smoother","None",jac->as_relax_type,&jac->as_relax_type,&flag);
1027:   PetscOptionsInt("-pc_hypre_ams_relax_times","Number of relaxation steps for AMS smoother","None",jac->as_relax_times,&jac->as_relax_times,&flag2);
1028:   PetscOptionsReal("-pc_hypre_ams_relax_weight","Relaxation weight for AMS smoother","None",jac->as_relax_weight,&jac->as_relax_weight,&flag3);
1029:   PetscOptionsReal("-pc_hypre_ams_omega","SSOR coefficient for AMS smoother","None",jac->as_omega,&jac->as_omega,&flag4);
1030:   if (flag || flag2 || flag3 || flag4) {
1031:     PetscStackCallStandard(HYPRE_AMSSetSmoothingOptions,(jac->hsolver,jac->as_relax_type,
1032:                                                                       jac->as_relax_times,
1033:                                                                       jac->as_relax_weight,
1034:                                                                       jac->as_omega));
1035:   }
1036:   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);
1037:   n = 5;
1038:   PetscOptionsIntArray("-pc_hypre_ams_amg_alpha_options","AMG options for vector Poisson","None",jac->as_amg_alpha_opts,&n,&flag2);
1039:   if (flag || flag2) {
1040:     PetscStackCallStandard(HYPRE_AMSSetAlphaAMGOptions,(jac->hsolver,jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
1041:                                                                      jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
1042:                                                                      jac->as_amg_alpha_opts[2],       /* AMG relax_type */
1043:                                                                      jac->as_amg_alpha_theta,
1044:                                                                      jac->as_amg_alpha_opts[3],       /* AMG interp_type */
1045:                                                                      jac->as_amg_alpha_opts[4]));     /* AMG Pmax */
1046:   }
1047:   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);
1048:   n = 5;
1049:   PetscOptionsIntArray("-pc_hypre_ams_amg_beta_options","AMG options for scalar Poisson solver","None",jac->as_amg_beta_opts,&n,&flag2);
1050:   if (flag || flag2) {
1051:     PetscStackCallStandard(HYPRE_AMSSetBetaAMGOptions,(jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
1052:                                                                     jac->as_amg_beta_opts[1],       /* AMG agg_levels */
1053:                                                                     jac->as_amg_beta_opts[2],       /* AMG relax_type */
1054:                                                                     jac->as_amg_beta_theta,
1055:                                                                     jac->as_amg_beta_opts[3],       /* AMG interp_type */
1056:                                                                     jac->as_amg_beta_opts[4]));     /* AMG Pmax */
1057:   }
1058:   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);
1059:   if (flag) { /* override HYPRE's default only if the options is used */
1060:     PetscStackCallStandard(HYPRE_AMSSetProjectionFrequency,(jac->hsolver,jac->ams_proj_freq));
1061:   }
1062:   PetscOptionsTail();
1063:   return(0);
1064: }

1066: static PetscErrorCode PCView_HYPRE_AMS(PC pc,PetscViewer viewer)
1067: {
1068:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1070:   PetscBool      iascii;

1073:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1074:   if (iascii) {
1075:     PetscViewerASCIIPrintf(viewer,"  HYPRE AMS preconditioning\n");
1076:     PetscViewerASCIIPrintf(viewer,"    subspace iterations per application %d\n",jac->as_max_iter);
1077:     PetscViewerASCIIPrintf(viewer,"    subspace cycle type %d\n",jac->ams_cycle_type);
1078:     PetscViewerASCIIPrintf(viewer,"    subspace iteration tolerance %g\n",jac->as_tol);
1079:     PetscViewerASCIIPrintf(viewer,"    smoother type %d\n",jac->as_relax_type);
1080:     PetscViewerASCIIPrintf(viewer,"    number of smoothing steps %d\n",jac->as_relax_times);
1081:     PetscViewerASCIIPrintf(viewer,"    smoother weight %g\n",jac->as_relax_weight);
1082:     PetscViewerASCIIPrintf(viewer,"    smoother omega %g\n",jac->as_omega);
1083:     if (jac->alpha_Poisson) {
1084:       PetscViewerASCIIPrintf(viewer,"    vector Poisson solver (passed in by user)\n");
1085:     } else {
1086:       PetscViewerASCIIPrintf(viewer,"    vector Poisson solver (computed) \n");
1087:     }
1088:     PetscViewerASCIIPrintf(viewer,"        boomerAMG coarsening type %d\n",jac->as_amg_alpha_opts[0]);
1089:     PetscViewerASCIIPrintf(viewer,"        boomerAMG levels of aggressive coarsening %d\n",jac->as_amg_alpha_opts[1]);
1090:     PetscViewerASCIIPrintf(viewer,"        boomerAMG relaxation type %d\n",jac->as_amg_alpha_opts[2]);
1091:     PetscViewerASCIIPrintf(viewer,"        boomerAMG interpolation type %d\n",jac->as_amg_alpha_opts[3]);
1092:     PetscViewerASCIIPrintf(viewer,"        boomerAMG max nonzero elements in interpolation rows %d\n",jac->as_amg_alpha_opts[4]);
1093:     PetscViewerASCIIPrintf(viewer,"        boomerAMG strength threshold %g\n",jac->as_amg_alpha_theta);
1094:     if (!jac->ams_beta_is_zero) {
1095:       if (jac->beta_Poisson) {
1096:         PetscViewerASCIIPrintf(viewer,"    scalar Poisson solver (passed in by user)\n");
1097:       } else {
1098:         PetscViewerASCIIPrintf(viewer,"    scalar Poisson solver (computed) \n");
1099:       }
1100:       PetscViewerASCIIPrintf(viewer,"        boomerAMG coarsening type %d\n",jac->as_amg_beta_opts[0]);
1101:       PetscViewerASCIIPrintf(viewer,"        boomerAMG levels of aggressive coarsening %d\n",jac->as_amg_beta_opts[1]);
1102:       PetscViewerASCIIPrintf(viewer,"        boomerAMG relaxation type %d\n",jac->as_amg_beta_opts[2]);
1103:       PetscViewerASCIIPrintf(viewer,"        boomerAMG interpolation type %d\n",jac->as_amg_beta_opts[3]);
1104:       PetscViewerASCIIPrintf(viewer,"        boomerAMG max nonzero elements in interpolation rows %d\n",jac->as_amg_beta_opts[4]);
1105:       PetscViewerASCIIPrintf(viewer,"        boomerAMG strength threshold %g\n",jac->as_amg_beta_theta);
1106:       if (jac->ams_beta_is_zero_part) {
1107:         PetscViewerASCIIPrintf(viewer,"        compatible subspace projection frequency %d (-1 HYPRE uses default)\n",jac->ams_proj_freq);
1108:       }
1109:     } else {
1110:       PetscViewerASCIIPrintf(viewer,"    scalar Poisson solver not used (zero-conductivity everywhere) \n");
1111:     }
1112:   }
1113:   return(0);
1114: }

1116: static PetscErrorCode PCSetFromOptions_HYPRE_ADS(PetscOptionItems *PetscOptionsObject,PC pc)
1117: {
1118:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1120:   PetscInt       n;
1121:   PetscBool      flag,flag2,flag3,flag4;

1124:   PetscOptionsHead(PetscOptionsObject,"HYPRE ADS Options");
1125:   PetscOptionsInt("-pc_hypre_ads_print_level","Debugging output level for ADS","None",jac->as_print,&jac->as_print,&flag);
1126:   if (flag) PetscStackCallStandard(HYPRE_ADSSetPrintLevel,(jac->hsolver,jac->as_print));
1127:   PetscOptionsInt("-pc_hypre_ads_max_iter","Maximum number of ADS multigrid iterations within PCApply","None",jac->as_max_iter,&jac->as_max_iter,&flag);
1128:   if (flag) PetscStackCallStandard(HYPRE_ADSSetMaxIter,(jac->hsolver,jac->as_max_iter));
1129:   PetscOptionsInt("-pc_hypre_ads_cycle_type","Cycle type for ADS multigrid","None",jac->ads_cycle_type,&jac->ads_cycle_type,&flag);
1130:   if (flag) PetscStackCallStandard(HYPRE_ADSSetCycleType,(jac->hsolver,jac->ads_cycle_type));
1131:   PetscOptionsReal("-pc_hypre_ads_tol","Error tolerance for ADS multigrid","None",jac->as_tol,&jac->as_tol,&flag);
1132:   if (flag) PetscStackCallStandard(HYPRE_ADSSetTol,(jac->hsolver,jac->as_tol));
1133:   PetscOptionsInt("-pc_hypre_ads_relax_type","Relaxation type for ADS smoother","None",jac->as_relax_type,&jac->as_relax_type,&flag);
1134:   PetscOptionsInt("-pc_hypre_ads_relax_times","Number of relaxation steps for ADS smoother","None",jac->as_relax_times,&jac->as_relax_times,&flag2);
1135:   PetscOptionsReal("-pc_hypre_ads_relax_weight","Relaxation weight for ADS smoother","None",jac->as_relax_weight,&jac->as_relax_weight,&flag3);
1136:   PetscOptionsReal("-pc_hypre_ads_omega","SSOR coefficient for ADS smoother","None",jac->as_omega,&jac->as_omega,&flag4);
1137:   if (flag || flag2 || flag3 || flag4) {
1138:     PetscStackCallStandard(HYPRE_ADSSetSmoothingOptions,(jac->hsolver,jac->as_relax_type,
1139:                                                                       jac->as_relax_times,
1140:                                                                       jac->as_relax_weight,
1141:                                                                       jac->as_omega));
1142:   }
1143:   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);
1144:   n = 5;
1145:   PetscOptionsIntArray("-pc_hypre_ads_ams_options","AMG options for AMS solver inside ADS","None",jac->as_amg_alpha_opts,&n,&flag2);
1146:   PetscOptionsInt("-pc_hypre_ads_ams_cycle_type","Cycle type for AMS solver inside ADS","None",jac->ams_cycle_type,&jac->ams_cycle_type,&flag3);
1147:   if (flag || flag2 || flag3) {
1148:     PetscStackCallStandard(HYPRE_ADSSetAMSOptions,(jac->hsolver,jac->ams_cycle_type,             /* AMS cycle type */
1149:                                                                 jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
1150:                                                                 jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
1151:                                                                 jac->as_amg_alpha_opts[2],       /* AMG relax_type */
1152:                                                                 jac->as_amg_alpha_theta,
1153:                                                                 jac->as_amg_alpha_opts[3],       /* AMG interp_type */
1154:                                                                 jac->as_amg_alpha_opts[4]));     /* AMG Pmax */
1155:   }
1156:   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);
1157:   n = 5;
1158:   PetscOptionsIntArray("-pc_hypre_ads_amg_options","AMG options for vector AMG solver inside ADS","None",jac->as_amg_beta_opts,&n,&flag2);
1159:   if (flag || flag2) {
1160:     PetscStackCallStandard(HYPRE_ADSSetAMGOptions,(jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
1161:                                                                 jac->as_amg_beta_opts[1],       /* AMG agg_levels */
1162:                                                                 jac->as_amg_beta_opts[2],       /* AMG relax_type */
1163:                                                                 jac->as_amg_beta_theta,
1164:                                                                 jac->as_amg_beta_opts[3],       /* AMG interp_type */
1165:                                                                 jac->as_amg_beta_opts[4]));     /* AMG Pmax */
1166:   }
1167:   PetscOptionsTail();
1168:   return(0);
1169: }

1171: static PetscErrorCode PCView_HYPRE_ADS(PC pc,PetscViewer viewer)
1172: {
1173:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1175:   PetscBool      iascii;

1178:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1179:   if (iascii) {
1180:     PetscViewerASCIIPrintf(viewer,"  HYPRE ADS preconditioning\n");
1181:     PetscViewerASCIIPrintf(viewer,"    subspace iterations per application %d\n",jac->as_max_iter);
1182:     PetscViewerASCIIPrintf(viewer,"    subspace cycle type %d\n",jac->ads_cycle_type);
1183:     PetscViewerASCIIPrintf(viewer,"    subspace iteration tolerance %g\n",jac->as_tol);
1184:     PetscViewerASCIIPrintf(viewer,"    smoother type %d\n",jac->as_relax_type);
1185:     PetscViewerASCIIPrintf(viewer,"    number of smoothing steps %d\n",jac->as_relax_times);
1186:     PetscViewerASCIIPrintf(viewer,"    smoother weight %g\n",jac->as_relax_weight);
1187:     PetscViewerASCIIPrintf(viewer,"    smoother omega %g\n",jac->as_omega);
1188:     PetscViewerASCIIPrintf(viewer,"    AMS solver using boomerAMG\n");
1189:     PetscViewerASCIIPrintf(viewer,"        subspace cycle type %d\n",jac->ams_cycle_type);
1190:     PetscViewerASCIIPrintf(viewer,"        coarsening type %d\n",jac->as_amg_alpha_opts[0]);
1191:     PetscViewerASCIIPrintf(viewer,"        levels of aggressive coarsening %d\n",jac->as_amg_alpha_opts[1]);
1192:     PetscViewerASCIIPrintf(viewer,"        relaxation type %d\n",jac->as_amg_alpha_opts[2]);
1193:     PetscViewerASCIIPrintf(viewer,"        interpolation type %d\n",jac->as_amg_alpha_opts[3]);
1194:     PetscViewerASCIIPrintf(viewer,"        max nonzero elements in interpolation rows %d\n",jac->as_amg_alpha_opts[4]);
1195:     PetscViewerASCIIPrintf(viewer,"        strength threshold %g\n",jac->as_amg_alpha_theta);
1196:     PetscViewerASCIIPrintf(viewer,"    vector Poisson solver using boomerAMG\n");
1197:     PetscViewerASCIIPrintf(viewer,"        coarsening type %d\n",jac->as_amg_beta_opts[0]);
1198:     PetscViewerASCIIPrintf(viewer,"        levels of aggressive coarsening %d\n",jac->as_amg_beta_opts[1]);
1199:     PetscViewerASCIIPrintf(viewer,"        relaxation type %d\n",jac->as_amg_beta_opts[2]);
1200:     PetscViewerASCIIPrintf(viewer,"        interpolation type %d\n",jac->as_amg_beta_opts[3]);
1201:     PetscViewerASCIIPrintf(viewer,"        max nonzero elements in interpolation rows %d\n",jac->as_amg_beta_opts[4]);
1202:     PetscViewerASCIIPrintf(viewer,"        strength threshold %g\n",jac->as_amg_beta_theta);
1203:   }
1204:   return(0);
1205: }

1207: static PetscErrorCode PCHYPRESetDiscreteGradient_HYPRE(PC pc, Mat G)
1208: {
1209:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1210:   PetscBool      ishypre;

1214:   PetscObjectTypeCompare((PetscObject)G,MATHYPRE,&ishypre);
1215:   if (ishypre) {
1216:     PetscObjectReference((PetscObject)G);
1217:     MatDestroy(&jac->G);
1218:     jac->G = G;
1219:   } else {
1220:     MatDestroy(&jac->G);
1221:     MatConvert(G,MATHYPRE,MAT_INITIAL_MATRIX,&jac->G);
1222:   }
1223:   return(0);
1224: }

1226: /*@
1227:  PCHYPRESetDiscreteGradient - Set discrete gradient matrix

1229:    Collective on PC

1231:    Input Parameters:
1232: +  pc - the preconditioning context
1233: -  G - the discrete gradient

1235:    Level: intermediate

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

1241: .seealso:
1242: @*/
1243: PetscErrorCode PCHYPRESetDiscreteGradient(PC pc, Mat G)
1244: {

1251:   PetscTryMethod(pc,"PCHYPRESetDiscreteGradient_C",(PC,Mat),(pc,G));
1252:   return(0);
1253: }

1255: static PetscErrorCode PCHYPRESetDiscreteCurl_HYPRE(PC pc, Mat C)
1256: {
1257:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1258:   PetscBool      ishypre;

1262:   PetscObjectTypeCompare((PetscObject)C,MATHYPRE,&ishypre);
1263:   if (ishypre) {
1264:     PetscObjectReference((PetscObject)C);
1265:     MatDestroy(&jac->C);
1266:     jac->C = C;
1267:   } else {
1268:     MatDestroy(&jac->C);
1269:     MatConvert(C,MATHYPRE,MAT_INITIAL_MATRIX,&jac->C);
1270:   }
1271:   return(0);
1272: }

1274: /*@
1275:  PCHYPRESetDiscreteCurl - Set discrete curl matrix

1277:    Collective on PC

1279:    Input Parameters:
1280: +  pc - the preconditioning context
1281: -  C - the discrete curl

1283:    Level: intermediate

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

1289: .seealso:
1290: @*/
1291: PetscErrorCode PCHYPRESetDiscreteCurl(PC pc, Mat C)
1292: {

1299:   PetscTryMethod(pc,"PCHYPRESetDiscreteCurl_C",(PC,Mat),(pc,C));
1300:   return(0);
1301: }

1303: static PetscErrorCode PCHYPRESetInterpolations_HYPRE(PC pc, PetscInt dim, Mat RT_PiFull, Mat RT_Pi[], Mat ND_PiFull, Mat ND_Pi[])
1304: {
1305:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1306:   PetscBool      ishypre;
1308:   PetscInt       i;

1311:   MatDestroy(&jac->RT_PiFull);
1312:   MatDestroy(&jac->ND_PiFull);
1313:   for (i=0;i<3;++i) {
1314:     MatDestroy(&jac->RT_Pi[i]);
1315:     MatDestroy(&jac->ND_Pi[i]);
1316:   }

1318:   jac->dim = dim;
1319:   if (RT_PiFull) {
1320:     PetscObjectTypeCompare((PetscObject)RT_PiFull,MATHYPRE,&ishypre);
1321:     if (ishypre) {
1322:       PetscObjectReference((PetscObject)RT_PiFull);
1323:       jac->RT_PiFull = RT_PiFull;
1324:     } else {
1325:       MatConvert(RT_PiFull,MATHYPRE,MAT_INITIAL_MATRIX,&jac->RT_PiFull);
1326:     }
1327:   }
1328:   if (RT_Pi) {
1329:     for (i=0;i<dim;++i) {
1330:       if (RT_Pi[i]) {
1331:         PetscObjectTypeCompare((PetscObject)RT_Pi[i],MATHYPRE,&ishypre);
1332:         if (ishypre) {
1333:           PetscObjectReference((PetscObject)RT_Pi[i]);
1334:           jac->RT_Pi[i] = RT_Pi[i];
1335:         } else {
1336:           MatConvert(RT_Pi[i],MATHYPRE,MAT_INITIAL_MATRIX,&jac->RT_Pi[i]);
1337:         }
1338:       }
1339:     }
1340:   }
1341:   if (ND_PiFull) {
1342:     PetscObjectTypeCompare((PetscObject)ND_PiFull,MATHYPRE,&ishypre);
1343:     if (ishypre) {
1344:       PetscObjectReference((PetscObject)ND_PiFull);
1345:       jac->ND_PiFull = ND_PiFull;
1346:     } else {
1347:       MatConvert(ND_PiFull,MATHYPRE,MAT_INITIAL_MATRIX,&jac->ND_PiFull);
1348:     }
1349:   }
1350:   if (ND_Pi) {
1351:     for (i=0;i<dim;++i) {
1352:       if (ND_Pi[i]) {
1353:         PetscObjectTypeCompare((PetscObject)ND_Pi[i],MATHYPRE,&ishypre);
1354:         if (ishypre) {
1355:           PetscObjectReference((PetscObject)ND_Pi[i]);
1356:           jac->ND_Pi[i] = ND_Pi[i];
1357:         } else {
1358:           MatConvert(ND_Pi[i],MATHYPRE,MAT_INITIAL_MATRIX,&jac->ND_Pi[i]);
1359:         }
1360:       }
1361:     }
1362:   }

1364:   return(0);
1365: }

1367: /*@
1368:  PCHYPRESetInterpolations - Set interpolation matrices for AMS/ADS preconditioner

1370:    Collective on PC

1372:    Input Parameters:
1373: +  pc - the preconditioning context
1374: -  dim - the dimension of the problem, only used in AMS
1375: -  RT_PiFull - Raviart-Thomas interpolation matrix
1376: -  RT_Pi - x/y/z component of Raviart-Thomas interpolation matrix
1377: -  ND_PiFull - Nedelec interpolation matrix
1378: -  ND_Pi - x/y/z component of Nedelec interpolation matrix

1380:    Notes:
1381:     For AMS, only Nedelec interpolation matrices are needed, the Raviart-Thomas interpolation matrices can be set to NULL.
1382:           For ADS, both type of interpolation matrices are needed.
1383:    Level: intermediate

1385: .seealso:
1386: @*/
1387: PetscErrorCode PCHYPRESetInterpolations(PC pc, PetscInt dim, Mat RT_PiFull, Mat RT_Pi[], Mat ND_PiFull, Mat ND_Pi[])
1388: {
1390:   PetscInt       i;

1394:   if (RT_PiFull) {
1397:   }
1398:   if (RT_Pi) {
1400:     for (i=0;i<dim;++i) {
1401:       if (RT_Pi[i]) {
1404:       }
1405:     }
1406:   }
1407:   if (ND_PiFull) {
1410:   }
1411:   if (ND_Pi) {
1413:     for (i=0;i<dim;++i) {
1414:       if (ND_Pi[i]) {
1417:       }
1418:     }
1419:   }
1420:   PetscTryMethod(pc,"PCHYPRESetInterpolations_C",(PC,PetscInt,Mat,Mat[],Mat,Mat[]),(pc,dim,RT_PiFull,RT_Pi,ND_PiFull,ND_Pi));
1421:   return(0);
1422: }

1424: static PetscErrorCode PCHYPRESetPoissonMatrix_HYPRE(PC pc, Mat A, PetscBool isalpha)
1425: {
1426:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1427:   PetscBool      ishypre;

1431:   PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);
1432:   if (ishypre) {
1433:     if (isalpha) {
1434:       PetscObjectReference((PetscObject)A);
1435:       MatDestroy(&jac->alpha_Poisson);
1436:       jac->alpha_Poisson = A;
1437:     } else {
1438:       if (A) {
1439:         PetscObjectReference((PetscObject)A);
1440:       } else {
1441:         jac->ams_beta_is_zero = PETSC_TRUE;
1442:       }
1443:       MatDestroy(&jac->beta_Poisson);
1444:       jac->beta_Poisson = A;
1445:     }
1446:   } else {
1447:     if (isalpha) {
1448:       MatDestroy(&jac->alpha_Poisson);
1449:       MatConvert(A,MATHYPRE,MAT_INITIAL_MATRIX,&jac->alpha_Poisson);
1450:     } else {
1451:       if (A) {
1452:         MatDestroy(&jac->beta_Poisson);
1453:         MatConvert(A,MATHYPRE,MAT_INITIAL_MATRIX,&jac->beta_Poisson);
1454:       } else {
1455:         MatDestroy(&jac->beta_Poisson);
1456:         jac->ams_beta_is_zero = PETSC_TRUE;
1457:       }
1458:     }
1459:   }
1460:   return(0);
1461: }

1463: /*@
1464:  PCHYPRESetAlphaPoissonMatrix - Set vector Poisson matrix

1466:    Collective on PC

1468:    Input Parameters:
1469: +  pc - the preconditioning context
1470: -  A - the matrix

1472:    Level: intermediate

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

1477: .seealso:
1478: @*/
1479: PetscErrorCode PCHYPRESetAlphaPoissonMatrix(PC pc, Mat A)
1480: {

1487:   PetscTryMethod(pc,"PCHYPRESetPoissonMatrix_C",(PC,Mat,PetscBool),(pc,A,PETSC_TRUE));
1488:   return(0);
1489: }

1491: /*@
1492:  PCHYPRESetBetaPoissonMatrix - Set Poisson matrix

1494:    Collective on PC

1496:    Input Parameters:
1497: +  pc - the preconditioning context
1498: -  A - the matrix

1500:    Level: intermediate

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

1506: .seealso:
1507: @*/
1508: PetscErrorCode PCHYPRESetBetaPoissonMatrix(PC pc, Mat A)
1509: {

1514:   if (A) {
1517:   }
1518:   PetscTryMethod(pc,"PCHYPRESetPoissonMatrix_C",(PC,Mat,PetscBool),(pc,A,PETSC_FALSE));
1519:   return(0);
1520: }

1522: static PetscErrorCode PCHYPRESetEdgeConstantVectors_HYPRE(PC pc,Vec ozz, Vec zoz, Vec zzo)
1523: {
1524:   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
1525:   PetscErrorCode     ierr;

1528:   /* throw away any vector if already set */
1529:   if (jac->constants[0]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[0]));
1530:   if (jac->constants[1]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[1]));
1531:   if (jac->constants[2]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->constants[2]));
1532:   jac->constants[0] = NULL;
1533:   jac->constants[1] = NULL;
1534:   jac->constants[2] = NULL;
1535:   VecHYPRE_IJVectorCreate(ozz,&jac->constants[0]);
1536:   VecHYPRE_IJVectorCopy(ozz,jac->constants[0]);
1537:   VecHYPRE_IJVectorCreate(zoz,&jac->constants[1]);
1538:   VecHYPRE_IJVectorCopy(zoz,jac->constants[1]);
1539:   jac->dim = 2;
1540:   if (zzo) {
1541:     VecHYPRE_IJVectorCreate(zzo,&jac->constants[2]);
1542:     VecHYPRE_IJVectorCopy(zzo,jac->constants[2]);
1543:     jac->dim++;
1544:   }
1545:   return(0);
1546: }

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

1551:    Collective on PC

1553:    Input Parameters:
1554: +  pc - the preconditioning context
1555: -  ozz - vector representing (1,0,0) (or (1,0) in 2D)
1556: -  zoz - vector representing (0,1,0) (or (0,1) in 2D)
1557: -  zzo - vector representing (0,0,1) (use NULL in 2D)

1559:    Level: intermediate

1561:    Notes:

1563: .seealso:
1564: @*/
1565: PetscErrorCode PCHYPRESetEdgeConstantVectors(PC pc, Vec ozz, Vec zoz, Vec zzo)
1566: {

1577:   PetscTryMethod(pc,"PCHYPRESetEdgeConstantVectors_C",(PC,Vec,Vec,Vec),(pc,ozz,zoz,zzo));
1578:   return(0);
1579: }

1581: static PetscErrorCode PCSetCoordinates_HYPRE(PC pc, PetscInt dim, PetscInt nloc, PetscReal *coords)
1582: {
1583:   PC_HYPRE        *jac = (PC_HYPRE*)pc->data;
1584:   Vec             tv;
1585:   PetscInt        i;
1586:   PetscErrorCode  ierr;

1589:   /* throw away any coordinate vector if already set */
1590:   if (jac->coords[0]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[0]));
1591:   if (jac->coords[1]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[1]));
1592:   if (jac->coords[2]) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->coords[2]));
1593:   jac->dim = dim;

1595:   /* compute IJ vector for coordinates */
1596:   VecCreate(PetscObjectComm((PetscObject)pc),&tv);
1597:   VecSetType(tv,VECSTANDARD);
1598:   VecSetSizes(tv,nloc,PETSC_DECIDE);
1599:   for (i=0;i<dim;i++) {
1600:     PetscScalar *array;
1601:     PetscInt    j;

1603:     VecHYPRE_IJVectorCreate(tv,&jac->coords[i]);
1604:     VecGetArray(tv,&array);
1605:     for (j=0;j<nloc;j++) {
1606:       array[j] = coords[j*dim+i];
1607:     }
1608:     PetscStackCallStandard(HYPRE_IJVectorSetValues,(jac->coords[i],nloc,NULL,array));
1609:     PetscStackCallStandard(HYPRE_IJVectorAssemble,(jac->coords[i]));
1610:     VecRestoreArray(tv,&array);
1611:   }
1612:   VecDestroy(&tv);
1613:   return(0);
1614: }

1616: /* ---------------------------------------------------------------------------------*/

1618: static PetscErrorCode  PCHYPREGetType_HYPRE(PC pc,const char *name[])
1619: {
1620:   PC_HYPRE *jac = (PC_HYPRE*)pc->data;

1623:   *name = jac->hypre_type;
1624:   return(0);
1625: }

1627: static PetscErrorCode  PCHYPRESetType_HYPRE(PC pc,const char name[])
1628: {
1629:   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
1631:   PetscBool      flag;

1634:   if (jac->hypre_type) {
1635:     PetscStrcmp(jac->hypre_type,name,&flag);
1636:     if (!flag) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
1637:     return(0);
1638:   } else {
1639:     PetscStrallocpy(name, &jac->hypre_type);
1640:   }

1642:   jac->maxiter         = PETSC_DEFAULT;
1643:   jac->tol             = PETSC_DEFAULT;
1644:   jac->printstatistics = PetscLogPrintInfo;

1646:   PetscStrcmp("pilut",jac->hypre_type,&flag);
1647:   if (flag) {
1648:     MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(jac->comm_hypre));
1649:     PetscStackCallStandard(HYPRE_ParCSRPilutCreate,(jac->comm_hypre,&jac->hsolver));
1650:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
1651:     pc->ops->view           = PCView_HYPRE_Pilut;
1652:     jac->destroy            = HYPRE_ParCSRPilutDestroy;
1653:     jac->setup              = HYPRE_ParCSRPilutSetup;
1654:     jac->solve              = HYPRE_ParCSRPilutSolve;
1655:     jac->factorrowsize      = PETSC_DEFAULT;
1656:     return(0);
1657:   }
1658:   PetscStrcmp("euclid",jac->hypre_type,&flag);
1659:   if (flag) {
1660:     MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(jac->comm_hypre));
1661:     PetscStackCallStandard(HYPRE_EuclidCreate,(jac->comm_hypre,&jac->hsolver));
1662:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
1663:     pc->ops->view           = PCView_HYPRE_Euclid;
1664:     jac->destroy            = HYPRE_EuclidDestroy;
1665:     jac->setup              = HYPRE_EuclidSetup;
1666:     jac->solve              = HYPRE_EuclidSolve;
1667:     jac->factorrowsize      = PETSC_DEFAULT;
1668:     jac->eu_level           = PETSC_DEFAULT; /* default */
1669:     return(0);
1670:   }
1671:   PetscStrcmp("parasails",jac->hypre_type,&flag);
1672:   if (flag) {
1673:     MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(jac->comm_hypre));
1674:     PetscStackCallStandard(HYPRE_ParaSailsCreate,(jac->comm_hypre,&jac->hsolver));
1675:     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
1676:     pc->ops->view           = PCView_HYPRE_ParaSails;
1677:     jac->destroy            = HYPRE_ParaSailsDestroy;
1678:     jac->setup              = HYPRE_ParaSailsSetup;
1679:     jac->solve              = HYPRE_ParaSailsSolve;
1680:     /* initialize */
1681:     jac->nlevels    = 1;
1682:     jac->threshhold = .1;
1683:     jac->filter     = .1;
1684:     jac->loadbal    = 0;
1685:     if (PetscLogPrintInfo) jac->logging = (int) PETSC_TRUE;
1686:     else jac->logging = (int) PETSC_FALSE;

1688:     jac->ruse = (int) PETSC_FALSE;
1689:     jac->symt = 0;
1690:     PetscStackCallStandard(HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels));
1691:     PetscStackCallStandard(HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter));
1692:     PetscStackCallStandard(HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal));
1693:     PetscStackCallStandard(HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging));
1694:     PetscStackCallStandard(HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse));
1695:     PetscStackCallStandard(HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
1696:     return(0);
1697:   }
1698:   PetscStrcmp("boomeramg",jac->hypre_type,&flag);
1699:   if (flag) {
1700:     HYPRE_BoomerAMGCreate(&jac->hsolver);
1701:     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
1702:     pc->ops->view            = PCView_HYPRE_BoomerAMG;
1703:     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
1704:     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
1705:     jac->destroy             = HYPRE_BoomerAMGDestroy;
1706:     jac->setup               = HYPRE_BoomerAMGSetup;
1707:     jac->solve               = HYPRE_BoomerAMGSolve;
1708:     jac->applyrichardson     = PETSC_FALSE;
1709:     /* these defaults match the hypre defaults */
1710:     jac->cycletype        = 1;
1711:     jac->maxlevels        = 25;
1712:     jac->maxiter          = 1;
1713:     jac->tol              = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */
1714:     jac->truncfactor      = 0.0;
1715:     jac->strongthreshold  = .25;
1716:     jac->maxrowsum        = .9;
1717:     jac->coarsentype      = 6;
1718:     jac->measuretype      = 0;
1719:     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = 1;
1720:     jac->smoothtype       = -1; /* Not set by default */
1721:     jac->smoothnumlevels  = 25;
1722:     jac->eu_level         = 0;
1723:     jac->eu_droptolerance = 0;
1724:     jac->eu_bj            = 0;
1725:     jac->relaxtype[0]     = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */
1726:     jac->relaxtype[2]     = 9; /*G.E. */
1727:     jac->relaxweight      = 1.0;
1728:     jac->outerrelaxweight = 1.0;
1729:     jac->relaxorder       = 1;
1730:     jac->interptype       = 0;
1731:     jac->agg_nl           = 0;
1732:     jac->pmax             = 0;
1733:     jac->truncfactor      = 0.0;
1734:     jac->agg_num_paths    = 1;

1736:     jac->nodal_coarsening      = 0;
1737:     jac->nodal_coarsening_diag = 0;
1738:     jac->vec_interp_variant    = 0;
1739:     jac->vec_interp_qmax       = 0;
1740:     jac->vec_interp_smooth     = PETSC_FALSE;
1741:     jac->interp_refine         = 0;
1742:     jac->nodal_relax           = PETSC_FALSE;
1743:     jac->nodal_relax_levels    = 1;
1744:     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
1745:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
1746:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
1747:     PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
1748:     PetscStackCallStandard(HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
1749:     PetscStackCallStandard(HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
1750:     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
1751:     PetscStackCallStandard(HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
1752:     PetscStackCallStandard(HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
1753:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
1754:     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
1755:     PetscStackCallStandard(HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
1756:     PetscStackCallStandard(HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
1757:     PetscStackCallStandard(HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
1758:     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, jac->relaxtype[0]));  /*defaults coarse to 9*/
1759:     PetscStackCallStandard(HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver, jac->gridsweeps[0])); /*defaults coarse to 1 */
1760:     return(0);
1761:   }
1762:   PetscStrcmp("ams",jac->hypre_type,&flag);
1763:   if (flag) {
1764:     HYPRE_AMSCreate(&jac->hsolver);
1765:     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_AMS;
1766:     pc->ops->view            = PCView_HYPRE_AMS;
1767:     jac->destroy             = HYPRE_AMSDestroy;
1768:     jac->setup               = HYPRE_AMSSetup;
1769:     jac->solve               = HYPRE_AMSSolve;
1770:     jac->coords[0]           = NULL;
1771:     jac->coords[1]           = NULL;
1772:     jac->coords[2]           = NULL;
1773:     /* solver parameters: these are borrowed from mfem package, and they are not the default values from HYPRE AMS */
1774:     jac->as_print           = 0;
1775:     jac->as_max_iter        = 1; /* used as a preconditioner */
1776:     jac->as_tol             = 0.; /* used as a preconditioner */
1777:     jac->ams_cycle_type     = 13;
1778:     /* Smoothing options */
1779:     jac->as_relax_type      = 2;
1780:     jac->as_relax_times     = 1;
1781:     jac->as_relax_weight    = 1.0;
1782:     jac->as_omega           = 1.0;
1783:     /* Vector valued Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
1784:     jac->as_amg_alpha_opts[0] = 10;
1785:     jac->as_amg_alpha_opts[1] = 1;
1786:     jac->as_amg_alpha_opts[2] = 6;
1787:     jac->as_amg_alpha_opts[3] = 6;
1788:     jac->as_amg_alpha_opts[4] = 4;
1789:     jac->as_amg_alpha_theta   = 0.25;
1790:     /* Scalar Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
1791:     jac->as_amg_beta_opts[0] = 10;
1792:     jac->as_amg_beta_opts[1] = 1;
1793:     jac->as_amg_beta_opts[2] = 6;
1794:     jac->as_amg_beta_opts[3] = 6;
1795:     jac->as_amg_beta_opts[4] = 4;
1796:     jac->as_amg_beta_theta   = 0.25;
1797:     PetscStackCallStandard(HYPRE_AMSSetPrintLevel,(jac->hsolver,jac->as_print));
1798:     PetscStackCallStandard(HYPRE_AMSSetMaxIter,(jac->hsolver,jac->as_max_iter));
1799:     PetscStackCallStandard(HYPRE_AMSSetCycleType,(jac->hsolver,jac->ams_cycle_type));
1800:     PetscStackCallStandard(HYPRE_AMSSetTol,(jac->hsolver,jac->as_tol));
1801:     PetscStackCallStandard(HYPRE_AMSSetSmoothingOptions,(jac->hsolver,jac->as_relax_type,
1802:                                                                       jac->as_relax_times,
1803:                                                                       jac->as_relax_weight,
1804:                                                                       jac->as_omega));
1805:     PetscStackCallStandard(HYPRE_AMSSetAlphaAMGOptions,(jac->hsolver,jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
1806:                                                                      jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
1807:                                                                      jac->as_amg_alpha_opts[2],       /* AMG relax_type */
1808:                                                                      jac->as_amg_alpha_theta,
1809:                                                                      jac->as_amg_alpha_opts[3],       /* AMG interp_type */
1810:                                                                      jac->as_amg_alpha_opts[4]));     /* AMG Pmax */
1811:     PetscStackCallStandard(HYPRE_AMSSetBetaAMGOptions,(jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
1812:                                                                     jac->as_amg_beta_opts[1],       /* AMG agg_levels */
1813:                                                                     jac->as_amg_beta_opts[2],       /* AMG relax_type */
1814:                                                                     jac->as_amg_beta_theta,
1815:                                                                     jac->as_amg_beta_opts[3],       /* AMG interp_type */
1816:                                                                     jac->as_amg_beta_opts[4]));     /* AMG Pmax */
1817:     /* Zero conductivity */
1818:     jac->ams_beta_is_zero      = PETSC_FALSE;
1819:     jac->ams_beta_is_zero_part = PETSC_FALSE;
1820:     return(0);
1821:   }
1822:   PetscStrcmp("ads",jac->hypre_type,&flag);
1823:   if (flag) {
1824:     HYPRE_ADSCreate(&jac->hsolver);
1825:     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_ADS;
1826:     pc->ops->view            = PCView_HYPRE_ADS;
1827:     jac->destroy             = HYPRE_ADSDestroy;
1828:     jac->setup               = HYPRE_ADSSetup;
1829:     jac->solve               = HYPRE_ADSSolve;
1830:     jac->coords[0]           = NULL;
1831:     jac->coords[1]           = NULL;
1832:     jac->coords[2]           = NULL;
1833:     /* solver parameters: these are borrowed from mfem package, and they are not the default values from HYPRE ADS */
1834:     jac->as_print           = 0;
1835:     jac->as_max_iter        = 1; /* used as a preconditioner */
1836:     jac->as_tol             = 0.; /* used as a preconditioner */
1837:     jac->ads_cycle_type     = 13;
1838:     /* Smoothing options */
1839:     jac->as_relax_type      = 2;
1840:     jac->as_relax_times     = 1;
1841:     jac->as_relax_weight    = 1.0;
1842:     jac->as_omega           = 1.0;
1843:     /* AMS solver parameters: cycle_type, coarsen type, agg_levels, relax_type, interp_type, Pmax */
1844:     jac->ams_cycle_type       = 14;
1845:     jac->as_amg_alpha_opts[0] = 10;
1846:     jac->as_amg_alpha_opts[1] = 1;
1847:     jac->as_amg_alpha_opts[2] = 6;
1848:     jac->as_amg_alpha_opts[3] = 6;
1849:     jac->as_amg_alpha_opts[4] = 4;
1850:     jac->as_amg_alpha_theta   = 0.25;
1851:     /* Vector Poisson AMG solver parameters: coarsen type, agg_levels, relax_type, interp_type, Pmax */
1852:     jac->as_amg_beta_opts[0] = 10;
1853:     jac->as_amg_beta_opts[1] = 1;
1854:     jac->as_amg_beta_opts[2] = 6;
1855:     jac->as_amg_beta_opts[3] = 6;
1856:     jac->as_amg_beta_opts[4] = 4;
1857:     jac->as_amg_beta_theta   = 0.25;
1858:     PetscStackCallStandard(HYPRE_ADSSetPrintLevel,(jac->hsolver,jac->as_print));
1859:     PetscStackCallStandard(HYPRE_ADSSetMaxIter,(jac->hsolver,jac->as_max_iter));
1860:     PetscStackCallStandard(HYPRE_ADSSetCycleType,(jac->hsolver,jac->ams_cycle_type));
1861:     PetscStackCallStandard(HYPRE_ADSSetTol,(jac->hsolver,jac->as_tol));
1862:     PetscStackCallStandard(HYPRE_ADSSetSmoothingOptions,(jac->hsolver,jac->as_relax_type,
1863:                                                                       jac->as_relax_times,
1864:                                                                       jac->as_relax_weight,
1865:                                                                       jac->as_omega));
1866:     PetscStackCallStandard(HYPRE_ADSSetAMSOptions,(jac->hsolver,jac->ams_cycle_type,             /* AMG coarsen type */
1867:                                                                 jac->as_amg_alpha_opts[0],       /* AMG coarsen type */
1868:                                                                 jac->as_amg_alpha_opts[1],       /* AMG agg_levels */
1869:                                                                 jac->as_amg_alpha_opts[2],       /* AMG relax_type */
1870:                                                                 jac->as_amg_alpha_theta,
1871:                                                                 jac->as_amg_alpha_opts[3],       /* AMG interp_type */
1872:                                                                 jac->as_amg_alpha_opts[4]));     /* AMG Pmax */
1873:     PetscStackCallStandard(HYPRE_ADSSetAMGOptions,(jac->hsolver,jac->as_amg_beta_opts[0],       /* AMG coarsen type */
1874:                                                                 jac->as_amg_beta_opts[1],       /* AMG agg_levels */
1875:                                                                 jac->as_amg_beta_opts[2],       /* AMG relax_type */
1876:                                                                 jac->as_amg_beta_theta,
1877:                                                                 jac->as_amg_beta_opts[3],       /* AMG interp_type */
1878:                                                                 jac->as_amg_beta_opts[4]));     /* AMG Pmax */
1879:     return(0);
1880:   }
1881:   PetscFree(jac->hypre_type);

1883:   jac->hypre_type = NULL;
1884:   SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are euclid, pilut, parasails, boomeramg, ams",name);
1885:   return(0);
1886: }

1888: /*
1889:     It only gets here if the HYPRE type has not been set before the call to
1890:    ...SetFromOptions() which actually is most of the time
1891: */
1892: static PetscErrorCode PCSetFromOptions_HYPRE(PetscOptionItems *PetscOptionsObject,PC pc)
1893: {
1895:   PetscInt       indx;
1896:   const char     *type[] = {"euclid","pilut","parasails","boomeramg","ams","ads"};
1897:   PetscBool      flg;

1900:   PetscOptionsHead(PetscOptionsObject,"HYPRE preconditioner options");
1901:   PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"boomeramg",&indx,&flg);
1902:   if (flg) {
1903:     PCHYPRESetType_HYPRE(pc,type[indx]);
1904:   } else {
1905:     PCHYPRESetType_HYPRE(pc,"boomeramg");
1906:   }
1907:   if (pc->ops->setfromoptions) {
1908:     pc->ops->setfromoptions(PetscOptionsObject,pc);
1909:   }
1910:   PetscOptionsTail();
1911:   return(0);
1912: }

1914: /*@C
1915:      PCHYPRESetType - Sets which hypre preconditioner you wish to use

1917:    Input Parameters:
1918: +     pc - the preconditioner context
1919: -     name - either  euclid, pilut, parasails, boomeramg, ams, ads

1921:    Options Database Keys:
1922:    -pc_hypre_type - One of euclid, pilut, parasails, boomeramg, ams, ads

1924:    Level: intermediate

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

1929: @*/
1930: PetscErrorCode  PCHYPRESetType(PC pc,const char name[])
1931: {

1937:   PetscTryMethod(pc,"PCHYPRESetType_C",(PC,const char[]),(pc,name));
1938:   return(0);
1939: }

1941: /*@C
1942:      PCHYPREGetType - Gets which hypre preconditioner you are using

1944:    Input Parameter:
1945: .     pc - the preconditioner context

1947:    Output Parameter:
1948: .     name - either  euclid, pilut, parasails, boomeramg, ams, ads

1950:    Level: intermediate

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

1955: @*/
1956: PetscErrorCode  PCHYPREGetType(PC pc,const char *name[])
1957: {

1963:   PetscTryMethod(pc,"PCHYPREGetType_C",(PC,const char*[]),(pc,name));
1964:   return(0);
1965: }

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

1970:    Options Database Keys:
1971: +   -pc_hypre_type - One of euclid, pilut, parasails, boomeramg, ams, ads
1972: -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
1973:           preconditioner

1975:    Level: intermediate

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

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

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

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

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

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

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

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

2012: M*/

2014: PETSC_EXTERN PetscErrorCode PCCreate_HYPRE(PC pc)
2015: {
2016:   PC_HYPRE       *jac;

2020:   PetscNewLog(pc,&jac);

2022:   pc->data                = jac;
2023:   pc->ops->reset          = PCReset_HYPRE;
2024:   pc->ops->destroy        = PCDestroy_HYPRE;
2025:   pc->ops->setfromoptions = PCSetFromOptions_HYPRE;
2026:   pc->ops->setup          = PCSetUp_HYPRE;
2027:   pc->ops->apply          = PCApply_HYPRE;
2028:   jac->comm_hypre         = MPI_COMM_NULL;
2029:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetType_C",PCHYPRESetType_HYPRE);
2030:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPREGetType_C",PCHYPREGetType_HYPRE);
2031:   PetscObjectComposeFunction((PetscObject)pc,"PCSetCoordinates_C",PCSetCoordinates_HYPRE);
2032:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteGradient_C",PCHYPRESetDiscreteGradient_HYPRE);
2033:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetDiscreteCurl_C",PCHYPRESetDiscreteCurl_HYPRE);
2034:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetInterpolations_C",PCHYPRESetInterpolations_HYPRE);
2035:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetEdgeConstantVectors_C",PCHYPRESetEdgeConstantVectors_HYPRE);
2036:   PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetPoissonMatrix_C",PCHYPRESetPoissonMatrix_HYPRE);
2037:   return(0);
2038: }

2040: /* ---------------------------------------------------------------------------------------------------------------------------------*/

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

2046:   /* keep copy of PFMG options used so may view them */
2047:   PetscInt its;
2048:   double   tol;
2049:   PetscInt relax_type;
2050:   PetscInt rap_type;
2051:   PetscInt num_pre_relax,num_post_relax;
2052:   PetscInt max_levels;
2053: } PC_PFMG;

2055: PetscErrorCode PCDestroy_PFMG(PC pc)
2056: {
2058:   PC_PFMG        *ex = (PC_PFMG*) pc->data;

2061:   if (ex->hsolver) PetscStackCallStandard(HYPRE_StructPFMGDestroy,(ex->hsolver));
2062:   MPI_Comm_free(&ex->hcomm);
2063:   PetscFree(pc->data);
2064:   return(0);
2065: }

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

2070: PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer)
2071: {
2073:   PetscBool      iascii;
2074:   PC_PFMG        *ex = (PC_PFMG*) pc->data;

2077:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
2078:   if (iascii) {
2079:     PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG preconditioning\n");
2080:     PetscViewerASCIIPrintf(viewer,"    max iterations %d\n",ex->its);
2081:     PetscViewerASCIIPrintf(viewer,"    tolerance %g\n",ex->tol);
2082:     PetscViewerASCIIPrintf(viewer,"    relax type %s\n",PFMGRelaxType[ex->relax_type]);
2083:     PetscViewerASCIIPrintf(viewer,"    RAP type %s\n",PFMGRAPType[ex->rap_type]);
2084:     PetscViewerASCIIPrintf(viewer,"    number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);
2085:     PetscViewerASCIIPrintf(viewer,"    max levels %d\n",ex->max_levels);
2086:   }
2087:   return(0);
2088: }

2090: PetscErrorCode PCSetFromOptions_PFMG(PetscOptionItems *PetscOptionsObject,PC pc)
2091: {
2093:   PC_PFMG        *ex = (PC_PFMG*) pc->data;
2094:   PetscBool      flg = PETSC_FALSE;

2097:   PetscOptionsHead(PetscOptionsObject,"PFMG options");
2098:   PetscOptionsBool("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",flg,&flg,NULL);
2099:   if (flg) {
2100:     int level=3;
2101:     PetscStackCallStandard(HYPRE_StructPFMGSetPrintLevel,(ex->hsolver,level));
2102:   }
2103:   PetscOptionsInt("-pc_pfmg_its","Number of iterations of PFMG to use as preconditioner","HYPRE_StructPFMGSetMaxIter",ex->its,&ex->its,NULL);
2104:   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(ex->hsolver,ex->its));
2105:   PetscOptionsInt("-pc_pfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_StructPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,NULL);
2106:   PetscStackCallStandard(HYPRE_StructPFMGSetNumPreRelax,(ex->hsolver,ex->num_pre_relax));
2107:   PetscOptionsInt("-pc_pfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_StructPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,NULL);
2108:   PetscStackCallStandard(HYPRE_StructPFMGSetNumPostRelax,(ex->hsolver,ex->num_post_relax));

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

2113:   PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,NULL);
2114:   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(ex->hsolver,ex->tol));
2115:   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);
2116:   PetscStackCallStandard(HYPRE_StructPFMGSetRelaxType,(ex->hsolver, ex->relax_type));
2117:   PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,ALEN(PFMGRAPType),PFMGRAPType[ex->rap_type],&ex->rap_type,NULL);
2118:   PetscStackCallStandard(HYPRE_StructPFMGSetRAPType,(ex->hsolver, ex->rap_type));
2119:   PetscOptionsTail();
2120:   return(0);
2121: }

2123: PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y)
2124: {
2125:   PetscErrorCode    ierr;
2126:   PC_PFMG           *ex = (PC_PFMG*) pc->data;
2127:   PetscScalar       *yy;
2128:   const PetscScalar *xx;
2129:   PetscInt          ilower[3],iupper[3];
2130:   Mat_HYPREStruct   *mx = (Mat_HYPREStruct*)(pc->pmat->data);

2133:   PetscCitationsRegister(hypreCitation,&cite);
2134:   DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);
2135:   iupper[0] += ilower[0] - 1;
2136:   iupper[1] += ilower[1] - 1;
2137:   iupper[2] += ilower[2] - 1;

2139:   /* copy x values over to hypre */
2140:   PetscStackCallStandard(HYPRE_StructVectorSetConstantValues,(mx->hb,0.0));
2141:   VecGetArrayRead(x,&xx);
2142:   PetscStackCallStandard(HYPRE_StructVectorSetBoxValues,(mx->hb,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,(PetscScalar*)xx));
2143:   VecRestoreArrayRead(x,&xx);
2144:   PetscStackCallStandard(HYPRE_StructVectorAssemble,(mx->hb));
2145:   PetscStackCallStandard(HYPRE_StructPFMGSolve,(ex->hsolver,mx->hmat,mx->hb,mx->hx));

2147:   /* copy solution values back to PETSc */
2148:   VecGetArray(y,&yy);
2149:   PetscStackCallStandard(HYPRE_StructVectorGetBoxValues,(mx->hx,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,yy));
2150:   VecRestoreArray(y,&yy);
2151:   return(0);
2152: }

2154: 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)
2155: {
2156:   PC_PFMG        *jac = (PC_PFMG*)pc->data;
2158:   PetscInt       oits;

2161:   PetscCitationsRegister(hypreCitation,&cite);
2162:   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(jac->hsolver,its*jac->its));
2163:   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(jac->hsolver,rtol));

2165:   PCApply_PFMG(pc,b,y);
2166:   PetscStackCallStandard(HYPRE_StructPFMGGetNumIterations,(jac->hsolver,(HYPRE_Int *)&oits));
2167:   *outits = oits;
2168:   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
2169:   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
2170:   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(jac->hsolver,jac->tol));
2171:   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(jac->hsolver,jac->its));
2172:   return(0);
2173: }


2176: PetscErrorCode PCSetUp_PFMG(PC pc)
2177: {
2178:   PetscErrorCode  ierr;
2179:   PC_PFMG         *ex = (PC_PFMG*) pc->data;
2180:   Mat_HYPREStruct *mx = (Mat_HYPREStruct*)(pc->pmat->data);
2181:   PetscBool       flg;

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

2187:   /* create the hypre solver object and set its information */
2188:   if (ex->hsolver) PetscStackCallStandard(HYPRE_StructPFMGDestroy,(ex->hsolver));
2189:   PetscStackCallStandard(HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
2190:   PetscStackCallStandard(HYPRE_StructPFMGSetup,(ex->hsolver,mx->hmat,mx->hb,mx->hx));
2191:   PetscStackCallStandard(HYPRE_StructPFMGSetZeroGuess,(ex->hsolver));
2192:   return(0);
2193: }

2195: /*MC
2196:      PCPFMG - the hypre PFMG multigrid solver

2198:    Level: advanced

2200:    Options Database:
2201: + -pc_pfmg_its <its> number of iterations of PFMG to use as preconditioner
2202: . -pc_pfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
2203: . -pc_pfmg_num_post_relax <steps> number of smoothing steps after coarse grid
2204: . -pc_pfmg_tol <tol> tolerance of PFMG
2205: . -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
2206: - -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin

2208:    Notes:
2209:     This is for CELL-centered descretizations

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

2214: .seealso:  PCMG, MATHYPRESTRUCT
2215: M*/

2217: PETSC_EXTERN PetscErrorCode PCCreate_PFMG(PC pc)
2218: {
2220:   PC_PFMG        *ex;

2223:   PetscNew(&ex); \
2224:   pc->data = ex;

2226:   ex->its            = 1;
2227:   ex->tol            = 1.e-8;
2228:   ex->relax_type     = 1;
2229:   ex->rap_type       = 0;
2230:   ex->num_pre_relax  = 1;
2231:   ex->num_post_relax = 1;
2232:   ex->max_levels     = 0;

2234:   pc->ops->setfromoptions  = PCSetFromOptions_PFMG;
2235:   pc->ops->view            = PCView_PFMG;
2236:   pc->ops->destroy         = PCDestroy_PFMG;
2237:   pc->ops->apply           = PCApply_PFMG;
2238:   pc->ops->applyrichardson = PCApplyRichardson_PFMG;
2239:   pc->ops->setup           = PCSetUp_PFMG;

2241:   MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ex->hcomm));
2242:   PetscStackCallStandard(HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
2243:   return(0);
2244: }

2246: /* ---------------------------------------------------------------------------------------------------------------------------------------------------*/

2248: /* we know we are working with a HYPRE_SStructMatrix */
2249: typedef struct {
2250:   MPI_Comm            hcomm;       /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */
2251:   HYPRE_SStructSolver ss_solver;

2253:   /* keep copy of SYSPFMG options used so may view them */
2254:   PetscInt its;
2255:   double   tol;
2256:   PetscInt relax_type;
2257:   PetscInt num_pre_relax,num_post_relax;
2258: } PC_SysPFMG;

2260: PetscErrorCode PCDestroy_SysPFMG(PC pc)
2261: {
2263:   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;

2266:   if (ex->ss_solver) PetscStackCallStandard(HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));
2267:   MPI_Comm_free(&ex->hcomm);
2268:   PetscFree(pc->data);
2269:   return(0);
2270: }

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

2274: PetscErrorCode PCView_SysPFMG(PC pc,PetscViewer viewer)
2275: {
2277:   PetscBool      iascii;
2278:   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;

2281:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
2282:   if (iascii) {
2283:     PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG preconditioning\n");
2284:     PetscViewerASCIIPrintf(viewer,"  max iterations %d\n",ex->its);
2285:     PetscViewerASCIIPrintf(viewer,"  tolerance %g\n",ex->tol);
2286:     PetscViewerASCIIPrintf(viewer,"  relax type %s\n",PFMGRelaxType[ex->relax_type]);
2287:     PetscViewerASCIIPrintf(viewer,"  number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);
2288:   }
2289:   return(0);
2290: }

2292: PetscErrorCode PCSetFromOptions_SysPFMG(PetscOptionItems *PetscOptionsObject,PC pc)
2293: {
2295:   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;
2296:   PetscBool      flg = PETSC_FALSE;

2299:   PetscOptionsHead(PetscOptionsObject,"SysPFMG options");
2300:   PetscOptionsBool("-pc_syspfmg_print_statistics","Print statistics","HYPRE_SStructSysPFMGSetPrintLevel",flg,&flg,NULL);
2301:   if (flg) {
2302:     int level=3;
2303:     PetscStackCallStandard(HYPRE_SStructSysPFMGSetPrintLevel,(ex->ss_solver,level));
2304:   }
2305:   PetscOptionsInt("-pc_syspfmg_its","Number of iterations of SysPFMG to use as preconditioner","HYPRE_SStructSysPFMGSetMaxIter",ex->its,&ex->its,NULL);
2306:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(ex->ss_solver,ex->its));
2307:   PetscOptionsInt("-pc_syspfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_SStructSysPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,NULL);
2308:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetNumPreRelax,(ex->ss_solver,ex->num_pre_relax));
2309:   PetscOptionsInt("-pc_syspfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_SStructSysPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,NULL);
2310:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetNumPostRelax,(ex->ss_solver,ex->num_post_relax));

2312:   PetscOptionsReal("-pc_syspfmg_tol","Tolerance of SysPFMG","HYPRE_SStructSysPFMGSetTol",ex->tol,&ex->tol,NULL);
2313:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(ex->ss_solver,ex->tol));
2314:   PetscOptionsEList("-pc_syspfmg_relax_type","Relax type for the up and down cycles","HYPRE_SStructSysPFMGSetRelaxType",SysPFMGRelaxType,ALEN(SysPFMGRelaxType),SysPFMGRelaxType[ex->relax_type],&ex->relax_type,NULL);
2315:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetRelaxType,(ex->ss_solver, ex->relax_type));
2316:   PetscOptionsTail();
2317:   return(0);
2318: }

2320: PetscErrorCode PCApply_SysPFMG(PC pc,Vec x,Vec y)
2321: {
2322:   PetscErrorCode    ierr;
2323:   PC_SysPFMG        *ex = (PC_SysPFMG*) pc->data;
2324:   PetscScalar       *yy;
2325:   const PetscScalar *xx;
2326:   PetscInt          ilower[3],iupper[3];
2327:   Mat_HYPRESStruct  *mx     = (Mat_HYPRESStruct*)(pc->pmat->data);
2328:   PetscInt          ordering= mx->dofs_order;
2329:   PetscInt          nvars   = mx->nvars;
2330:   PetscInt          part    = 0;
2331:   PetscInt          size;
2332:   PetscInt          i;

2335:   PetscCitationsRegister(hypreCitation,&cite);
2336:   DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);
2337:   iupper[0] += ilower[0] - 1;
2338:   iupper[1] += ilower[1] - 1;
2339:   iupper[2] += ilower[2] - 1;

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

2344:   /* copy x values over to hypre for variable ordering */
2345:   if (ordering) {
2346:     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
2347:     VecGetArrayRead(x,&xx);
2348:     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,(PetscScalar*)xx+(size*i)));
2349:     VecRestoreArrayRead(x,&xx);
2350:     PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
2351:     PetscStackCallStandard(HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
2352:     PetscStackCallStandard(HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));

2354:     /* copy solution values back to PETSc */
2355:     VecGetArray(y,&yy);
2356:     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,yy+(size*i)));
2357:     VecRestoreArray(y,&yy);
2358:   } else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
2359:     PetscScalar *z;
2360:     PetscInt    j, k;

2362:     PetscMalloc1(nvars*size,&z);
2363:     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
2364:     VecGetArrayRead(x,&xx);

2366:     /* transform nodal to hypre's variable ordering for sys_pfmg */
2367:     for (i= 0; i< size; i++) {
2368:       k= i*nvars;
2369:       for (j= 0; j< nvars; j++) z[j*size+i]= xx[k+j];
2370:     }
2371:     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,z+(size*i)));
2372:     VecRestoreArrayRead(x,&xx);
2373:     PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
2374:     PetscStackCallStandard(HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));

2376:     /* copy solution values back to PETSc */
2377:     VecGetArray(y,&yy);
2378:     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,z+(size*i)));
2379:     /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
2380:     for (i= 0; i< size; i++) {
2381:       k= i*nvars;
2382:       for (j= 0; j< nvars; j++) yy[k+j]= z[j*size+i];
2383:     }
2384:     VecRestoreArray(y,&yy);
2385:     PetscFree(z);
2386:   }
2387:   return(0);
2388: }

2390: 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)
2391: {
2392:   PC_SysPFMG     *jac = (PC_SysPFMG*)pc->data;
2394:   PetscInt       oits;

2397:   PetscCitationsRegister(hypreCitation,&cite);
2398:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,its*jac->its));
2399:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,rtol));
2400:   PCApply_SysPFMG(pc,b,y);
2401:   PetscStackCallStandard(HYPRE_SStructSysPFMGGetNumIterations,(jac->ss_solver,(HYPRE_Int *)&oits));
2402:   *outits = oits;
2403:   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
2404:   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
2405:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,jac->tol));
2406:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,jac->its));
2407:   return(0);
2408: }

2410: PetscErrorCode PCSetUp_SysPFMG(PC pc)
2411: {
2412:   PetscErrorCode   ierr;
2413:   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
2414:   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct*)(pc->pmat->data);
2415:   PetscBool        flg;

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

2421:   /* create the hypre sstruct solver object and set its information */
2422:   if (ex->ss_solver) PetscStackCallStandard(HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));
2423:   PetscStackCallStandard(HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
2424:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetZeroGuess,(ex->ss_solver));
2425:   PetscStackCallStandard(HYPRE_SStructSysPFMGSetup,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
2426:   return(0);
2427: }

2429: /*MC
2430:      PCSysPFMG - the hypre SysPFMG multigrid solver

2432:    Level: advanced

2434:    Options Database:
2435: + -pc_syspfmg_its <its> number of iterations of SysPFMG to use as preconditioner
2436: . -pc_syspfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
2437: . -pc_syspfmg_num_post_relax <steps> number of smoothing steps after coarse grid
2438: . -pc_syspfmg_tol <tol> tolerance of SysPFMG
2439: . -pc_syspfmg_relax_type -relaxation type for the up and down cycles, one of Weighted-Jacobi,Red/Black-Gauss-Seidel

2441:    Notes:
2442:     This is for CELL-centered descretizations

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

2448: .seealso:  PCMG, MATHYPRESSTRUCT
2449: M*/

2451: PETSC_EXTERN PetscErrorCode PCCreate_SysPFMG(PC pc)
2452: {
2454:   PC_SysPFMG     *ex;

2457:   PetscNew(&ex); \
2458:   pc->data = ex;

2460:   ex->its            = 1;
2461:   ex->tol            = 1.e-8;
2462:   ex->relax_type     = 1;
2463:   ex->num_pre_relax  = 1;
2464:   ex->num_post_relax = 1;

2466:   pc->ops->setfromoptions  = PCSetFromOptions_SysPFMG;
2467:   pc->ops->view            = PCView_SysPFMG;
2468:   pc->ops->destroy         = PCDestroy_SysPFMG;
2469:   pc->ops->apply           = PCApply_SysPFMG;
2470:   pc->ops->applyrichardson = PCApplyRichardson_SysPFMG;
2471:   pc->ops->setup           = PCSetUp_SysPFMG;

2473:   MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ex->hcomm));
2474:   PetscStackCallStandard(HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
2475:   return(0);
2476: }