xref: /petsc/src/ksp/pc/impls/hypre/hypre.c (revision 8cc058d9cd56c1ccb3be12a47760ddfc446aaffc)
116d9e3a6SLisandro Dalcin 
216d9e3a6SLisandro Dalcin /*
316d9e3a6SLisandro Dalcin    Provides an interface to the LLNL package hypre
416d9e3a6SLisandro Dalcin */
50f1074feSSatish Balay 
60f1074feSSatish Balay /* Must use hypre 2.0.0 or more recent. */
70f1074feSSatish Balay 
8b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>          /*I "petscpc.h" I*/
9c6db04a5SJed Brown #include <../src/dm/impls/da/hypre/mhyp.h>
1016d9e3a6SLisandro Dalcin 
1116d9e3a6SLisandro Dalcin /*
1216d9e3a6SLisandro Dalcin    Private context (data structure) for the  preconditioner.
1316d9e3a6SLisandro Dalcin */
1416d9e3a6SLisandro Dalcin typedef struct {
1516d9e3a6SLisandro Dalcin   HYPRE_Solver   hsolver;
1616d9e3a6SLisandro Dalcin   HYPRE_IJMatrix ij;
1716d9e3a6SLisandro Dalcin   HYPRE_IJVector b,x;
1816d9e3a6SLisandro Dalcin 
194ddd07fcSJed Brown   HYPRE_Int (*destroy)(HYPRE_Solver);
204ddd07fcSJed Brown   HYPRE_Int (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
214ddd07fcSJed Brown   HYPRE_Int (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
2216d9e3a6SLisandro Dalcin 
2316d9e3a6SLisandro Dalcin   MPI_Comm comm_hypre;
2416d9e3a6SLisandro Dalcin   char     *hypre_type;
2516d9e3a6SLisandro Dalcin 
2616d9e3a6SLisandro Dalcin   /* options for Pilut and BoomerAMG*/
274ddd07fcSJed Brown   PetscInt maxiter;
2816d9e3a6SLisandro Dalcin   double   tol;
2916d9e3a6SLisandro Dalcin 
3016d9e3a6SLisandro Dalcin   /* options for Pilut */
314ddd07fcSJed Brown   PetscInt factorrowsize;
3216d9e3a6SLisandro Dalcin 
3316d9e3a6SLisandro Dalcin   /* options for ParaSails */
344ddd07fcSJed Brown   PetscInt nlevels;
3516d9e3a6SLisandro Dalcin   double   threshhold;
3616d9e3a6SLisandro Dalcin   double   filter;
374ddd07fcSJed Brown   PetscInt sym;
3816d9e3a6SLisandro Dalcin   double   loadbal;
394ddd07fcSJed Brown   PetscInt logging;
404ddd07fcSJed Brown   PetscInt ruse;
414ddd07fcSJed Brown   PetscInt symt;
4216d9e3a6SLisandro Dalcin 
4316d9e3a6SLisandro Dalcin   /* options for Euclid */
44ace3abfcSBarry Smith   PetscBool bjilu;
454ddd07fcSJed Brown   PetscInt  levels;
4616d9e3a6SLisandro Dalcin 
4716d9e3a6SLisandro Dalcin   /* options for Euclid and BoomerAMG */
48ace3abfcSBarry Smith   PetscBool printstatistics;
4916d9e3a6SLisandro Dalcin 
5016d9e3a6SLisandro Dalcin   /* options for BoomerAMG */
514ddd07fcSJed Brown   PetscInt  cycletype;
524ddd07fcSJed Brown   PetscInt  maxlevels;
5316d9e3a6SLisandro Dalcin   double    strongthreshold;
5416d9e3a6SLisandro Dalcin   double    maxrowsum;
554ddd07fcSJed Brown   PetscInt  gridsweeps[3];
564ddd07fcSJed Brown   PetscInt  coarsentype;
574ddd07fcSJed Brown   PetscInt  measuretype;
584ddd07fcSJed Brown   PetscInt  relaxtype[3];
5916d9e3a6SLisandro Dalcin   double    relaxweight;
6016d9e3a6SLisandro Dalcin   double    outerrelaxweight;
614ddd07fcSJed Brown   PetscInt  relaxorder;
6216d9e3a6SLisandro Dalcin   double    truncfactor;
63ace3abfcSBarry Smith   PetscBool applyrichardson;
644ddd07fcSJed Brown   PetscInt  pmax;
654ddd07fcSJed Brown   PetscInt  interptype;
664ddd07fcSJed Brown   PetscInt  agg_nl;
674ddd07fcSJed Brown   PetscInt  agg_num_paths;
684ddd07fcSJed Brown   PetscInt  nodal_coarsen;
69ace3abfcSBarry Smith   PetscBool nodal_relax;
704ddd07fcSJed Brown   PetscInt  nodal_relax_levels;
7116d9e3a6SLisandro Dalcin } PC_HYPRE;
7216d9e3a6SLisandro Dalcin 
73d2128fa2SBarry Smith #undef __FUNCT__
74d2128fa2SBarry Smith #define __FUNCT__ "PCHYPREGetSolver"
75d2128fa2SBarry Smith PetscErrorCode PCHYPREGetSolver(PC pc,HYPRE_Solver *hsolver)
76d2128fa2SBarry Smith {
77d2128fa2SBarry Smith   PC_HYPRE *jac = (PC_HYPRE*)pc->data;
78d2128fa2SBarry Smith 
79d2128fa2SBarry Smith   PetscFunctionBegin;
80d2128fa2SBarry Smith   *hsolver = jac->hsolver;
81d2128fa2SBarry Smith   PetscFunctionReturn(0);
82d2128fa2SBarry Smith }
8316d9e3a6SLisandro Dalcin 
8416d9e3a6SLisandro Dalcin #undef __FUNCT__
8516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetUp_HYPRE"
8616d9e3a6SLisandro Dalcin static PetscErrorCode PCSetUp_HYPRE(PC pc)
8716d9e3a6SLisandro Dalcin {
8816d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
8916d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
9016d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
9116d9e3a6SLisandro Dalcin   HYPRE_ParVector    bv,xv;
9216d9e3a6SLisandro Dalcin   PetscInt           bs;
9316d9e3a6SLisandro Dalcin 
9416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9516d9e3a6SLisandro Dalcin   if (!jac->hypre_type) {
9602a17cd4SBarry Smith     ierr = PCHYPRESetType(pc,"boomeramg");CHKERRQ(ierr);
9716d9e3a6SLisandro Dalcin   }
985f5c5b43SBarry Smith 
995f5c5b43SBarry Smith   if (pc->setupcalled) {
1005f5c5b43SBarry Smith     /* always destroy the old matrix and create a new memory;
1015f5c5b43SBarry Smith        hope this does not churn the memory too much. The problem
1025f5c5b43SBarry Smith        is I do not know if it is possible to put the matrix back to
1035f5c5b43SBarry Smith        its initial state so that we can directly copy the values
1045f5c5b43SBarry Smith        the second time through. */
105fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->ij));
1065f5c5b43SBarry Smith     jac->ij = 0;
10716d9e3a6SLisandro Dalcin   }
1085f5c5b43SBarry Smith 
10916d9e3a6SLisandro Dalcin   if (!jac->ij) { /* create the matrix the first time through */
11016d9e3a6SLisandro Dalcin     ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr);
11116d9e3a6SLisandro Dalcin   }
11216d9e3a6SLisandro Dalcin   if (!jac->b) { /* create the vectors the first time through */
11316d9e3a6SLisandro Dalcin     Vec x,b;
11416d9e3a6SLisandro Dalcin     ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr);
11516d9e3a6SLisandro Dalcin     ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr);
11616d9e3a6SLisandro Dalcin     ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr);
1176bf464f9SBarry Smith     ierr = VecDestroy(&x);CHKERRQ(ierr);
1186bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
11916d9e3a6SLisandro Dalcin   }
1205f5c5b43SBarry Smith 
12116d9e3a6SLisandro Dalcin   /* special case for BoomerAMG */
12216d9e3a6SLisandro Dalcin   if (jac->setup == HYPRE_BoomerAMGSetup) {
12316d9e3a6SLisandro Dalcin     ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr);
1242fa5cd67SKarl Rupp     if (bs > 1) PetscStackCallStandard(HYPRE_BoomerAMGSetNumFunctions,(jac->hsolver,bs));
12516d9e3a6SLisandro Dalcin   };
12616d9e3a6SLisandro Dalcin   ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr);
127fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
128fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&bv));
129fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&xv));
130fd3f9acdSBarry Smith   PetscStackCall("HYPRE_SetupXXX",ierr = (*jac->setup)(jac->hsolver,hmat,bv,xv);CHKERRQ(ierr););
13116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
13216d9e3a6SLisandro Dalcin }
13316d9e3a6SLisandro Dalcin 
13416d9e3a6SLisandro Dalcin /*
13516d9e3a6SLisandro Dalcin     Replaces the address where the HYPRE vector points to its data with the address of
13616d9e3a6SLisandro Dalcin   PETSc's data. Saves the old address so it can be reset when we are finished with it.
13716d9e3a6SLisandro Dalcin   Allows use to get the data into a HYPRE vector without the cost of memcopies
13816d9e3a6SLisandro Dalcin */
13916d9e3a6SLisandro Dalcin #define HYPREReplacePointer(b,newvalue,savedvalue) { \
14016d9e3a6SLisandro Dalcin     hypre_ParVector *par_vector   = (hypre_ParVector*)hypre_IJVectorObject(((hypre_IJVector*)b)); \
14116d9e3a6SLisandro Dalcin     hypre_Vector    *local_vector = hypre_ParVectorLocalVector(par_vector); \
14216d9e3a6SLisandro Dalcin     savedvalue         = local_vector->data; \
1430ad7597dSKarl Rupp     local_vector->data = newvalue;          \
1440ad7597dSKarl Rupp }
14516d9e3a6SLisandro Dalcin 
14616d9e3a6SLisandro Dalcin #undef __FUNCT__
14716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApply_HYPRE"
14816d9e3a6SLisandro Dalcin static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
14916d9e3a6SLisandro Dalcin {
15016d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
15116d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
15216d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
15316d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
15416d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
15516d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
1564ddd07fcSJed Brown   PetscInt           hierr;
15716d9e3a6SLisandro Dalcin 
15816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
15916d9e3a6SLisandro Dalcin   if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);}
16016d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
16116d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
16216d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
16316d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
16416d9e3a6SLisandro Dalcin 
165fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
166fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&jbv));
167fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&jxv));
168fd3f9acdSBarry Smith   PetscStackCall("Hypre solve",hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv);
16965e19b50SBarry Smith                                if (hierr && hierr != HYPRE_ERROR_CONV) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
170fd3f9acdSBarry Smith                                if (hierr) hypre__global_error = 0;);
17116d9e3a6SLisandro Dalcin 
17216d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
17316d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
17416d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
17516d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
17616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
17716d9e3a6SLisandro Dalcin }
17816d9e3a6SLisandro Dalcin 
17916d9e3a6SLisandro Dalcin #undef __FUNCT__
18016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCDestroy_HYPRE"
18116d9e3a6SLisandro Dalcin static PetscErrorCode PCDestroy_HYPRE(PC pc)
18216d9e3a6SLisandro Dalcin {
18316d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
18416d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
18516d9e3a6SLisandro Dalcin 
18616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
187fd3f9acdSBarry Smith   if (jac->ij) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->ij));
188fd3f9acdSBarry Smith   if (jac->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->b));
189fd3f9acdSBarry Smith   if (jac->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->x));
190fd3f9acdSBarry Smith   if (jac->destroy) PetscStackCall("HYPRE_DistroyXXX",ierr = (*jac->destroy)(jac->hsolver);CHKERRQ(ierr););
191503cfb0cSBarry Smith   ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr);
19216d9e3a6SLisandro Dalcin   if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);}
193c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
19416d9e3a6SLisandro Dalcin 
19516d9e3a6SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr);
19600de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetType_C","",NULL);CHKERRQ(ierr);
19700de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHYPREGetType_C","",NULL);CHKERRQ(ierr);
19816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
19916d9e3a6SLisandro Dalcin }
20016d9e3a6SLisandro Dalcin 
20116d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
20216d9e3a6SLisandro Dalcin #undef __FUNCT__
20316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut"
20416d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc)
20516d9e3a6SLisandro Dalcin {
20616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
20716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
208ace3abfcSBarry Smith   PetscBool      flag;
20916d9e3a6SLisandro Dalcin 
21016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
21116d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr);
21216d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr);
213fd3f9acdSBarry Smith   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetMaxIter,(jac->hsolver,jac->maxiter));
21416d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr);
215fd3f9acdSBarry Smith   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetDropTolerance,(jac->hsolver,jac->tol));
21616d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr);
217fd3f9acdSBarry Smith   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetFactorRowSize,(jac->hsolver,jac->factorrowsize));
21816d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
21916d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
22016d9e3a6SLisandro Dalcin }
22116d9e3a6SLisandro Dalcin 
22216d9e3a6SLisandro Dalcin #undef __FUNCT__
22316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Pilut"
22416d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
22516d9e3a6SLisandro Dalcin {
22616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
22716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
228ace3abfcSBarry Smith   PetscBool      iascii;
22916d9e3a6SLisandro Dalcin 
23016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
231251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
23216d9e3a6SLisandro Dalcin   if (iascii) {
23316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");CHKERRQ(ierr);
23416d9e3a6SLisandro Dalcin     if (jac->maxiter != PETSC_DEFAULT) {
23516d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr);
23616d9e3a6SLisandro Dalcin     } else {
23716d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr);
23816d9e3a6SLisandro Dalcin     }
23916d9e3a6SLisandro Dalcin     if (jac->tol != PETSC_DEFAULT) {
24016d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr);
24116d9e3a6SLisandro Dalcin     } else {
24216d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr);
24316d9e3a6SLisandro Dalcin     }
24416d9e3a6SLisandro Dalcin     if (jac->factorrowsize != PETSC_DEFAULT) {
24516d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr);
24616d9e3a6SLisandro Dalcin     } else {
24716d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");CHKERRQ(ierr);
24816d9e3a6SLisandro Dalcin     }
24916d9e3a6SLisandro Dalcin   }
25016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
25116d9e3a6SLisandro Dalcin }
25216d9e3a6SLisandro Dalcin 
25316d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
25416d9e3a6SLisandro Dalcin #undef __FUNCT__
25516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid"
25616d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc)
25716d9e3a6SLisandro Dalcin {
25816d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
25916d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
260ace3abfcSBarry Smith   PetscBool      flag;
261390e7148SBarry Smith   char           *args[8],levels[16];
262390e7148SBarry Smith   PetscInt       cnt = 0;
26316d9e3a6SLisandro Dalcin 
26416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
26516d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr);
26616d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr);
26716d9e3a6SLisandro Dalcin   if (flag) {
268ce94432eSBarry Smith     if (jac->levels < 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels);
2698caf3d72SBarry Smith     ierr        = PetscSNPrintf(levels,sizeof(levels),"%D",jac->levels);CHKERRQ(ierr);
270390e7148SBarry Smith     args[cnt++] = (char*)"-level"; args[cnt++] = levels;
27116d9e3a6SLisandro Dalcin   }
2720298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,NULL);CHKERRQ(ierr);
27316d9e3a6SLisandro Dalcin   if (jac->bjilu) {
274390e7148SBarry Smith     args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1";
27516d9e3a6SLisandro Dalcin   }
27616d9e3a6SLisandro Dalcin 
2770298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,NULL);CHKERRQ(ierr);
27816d9e3a6SLisandro Dalcin   if (jac->printstatistics) {
279390e7148SBarry Smith     args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1";
280390e7148SBarry Smith     args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1";
28116d9e3a6SLisandro Dalcin   }
28216d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
283fd3f9acdSBarry Smith   if (cnt) PetscStackCallStandard(HYPRE_EuclidSetParams,(jac->hsolver,cnt,args));
28416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
28516d9e3a6SLisandro Dalcin }
28616d9e3a6SLisandro Dalcin 
28716d9e3a6SLisandro Dalcin #undef __FUNCT__
28816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Euclid"
28916d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
29016d9e3a6SLisandro Dalcin {
29116d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
29216d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
293ace3abfcSBarry Smith   PetscBool      iascii;
29416d9e3a6SLisandro Dalcin 
29516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
296251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
29716d9e3a6SLisandro Dalcin   if (iascii) {
29816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n");CHKERRQ(ierr);
29916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr);
30016d9e3a6SLisandro Dalcin     if (jac->bjilu) {
30116d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr);
30216d9e3a6SLisandro Dalcin     }
30316d9e3a6SLisandro Dalcin   }
30416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
30516d9e3a6SLisandro Dalcin }
30616d9e3a6SLisandro Dalcin 
30716d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
30816d9e3a6SLisandro Dalcin 
30916d9e3a6SLisandro Dalcin #undef __FUNCT__
31016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG"
31116d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
31216d9e3a6SLisandro Dalcin {
31316d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
31416d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
31516d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
31616d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
31716d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
31816d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
3194ddd07fcSJed Brown   PetscInt           hierr;
32016d9e3a6SLisandro Dalcin 
32116d9e3a6SLisandro Dalcin   PetscFunctionBegin;
32216d9e3a6SLisandro Dalcin   ierr = VecSet(x,0.0);CHKERRQ(ierr);
32316d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
32416d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
32516d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
32616d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
32716d9e3a6SLisandro Dalcin 
328fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
329fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&jbv));
330fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&jxv));
33116d9e3a6SLisandro Dalcin 
33216d9e3a6SLisandro Dalcin   hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv);
33316d9e3a6SLisandro Dalcin   /* error code of 1 in BoomerAMG merely means convergence not achieved */
334e32f2f54SBarry Smith   if (hierr && (hierr != 1)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
33516d9e3a6SLisandro Dalcin   if (hierr) hypre__global_error = 0;
33616d9e3a6SLisandro Dalcin 
33716d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
33816d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
33916d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
34016d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
34116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
34216d9e3a6SLisandro Dalcin }
34316d9e3a6SLisandro Dalcin 
344a669f990SJed Brown /* static array length */
345a669f990SJed Brown #define ALEN(a) (sizeof(a)/sizeof((a)[0]))
346a669f990SJed Brown 
34716d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
3480f1074feSSatish Balay static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"};
34916d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
35065de4495SJed Brown /* The following corresponds to HYPRE_BoomerAMGSetRelaxType which has many missing numbers in the enum */
35165de4495SJed Brown static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","seqboundary-Gauss-Seidel","SOR/Jacobi","backward-SOR/Jacobi",
35265de4495SJed Brown                                                   "" /* [5] hybrid chaotic Gauss-Seidel (works only with OpenMP) */,"symmetric-SOR/Jacobi",
35365de4495SJed Brown                                                   "" /* 7 */,"l1scaled-SOR/Jacobi","Gaussian-elimination",
35465de4495SJed Brown                                                   "" /* 10 */, "" /* 11 */, "" /* 12 */, "" /* 13 */, "" /* 14 */,
35565de4495SJed Brown                                                   "CG" /* non-stationary */,"Chebyshev","FCF-Jacobi","l1scaled-Jacobi"};
3560f1074feSSatish Balay static const char *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i",
3570f1074feSSatish Balay                                                   "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"};
35816d9e3a6SLisandro Dalcin #undef __FUNCT__
35916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG"
36016d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
36116d9e3a6SLisandro Dalcin {
36216d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
36316d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
3644ddd07fcSJed Brown   PetscInt       n,indx,level;
365ace3abfcSBarry Smith   PetscBool      flg, tmp_truth;
36616d9e3a6SLisandro Dalcin   double         tmpdbl, twodbl[2];
36716d9e3a6SLisandro Dalcin 
36816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
36916d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr);
3704336a9eeSBarry Smith   ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType+1,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr);
37116d9e3a6SLisandro Dalcin   if (flg) {
3724336a9eeSBarry Smith     jac->cycletype = indx+1;
373fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
37416d9e3a6SLisandro Dalcin   }
37516d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr);
37616d9e3a6SLisandro Dalcin   if (flg) {
377ce94432eSBarry Smith     if (jac->maxlevels < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
378fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
37916d9e3a6SLisandro Dalcin   }
38016d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr);
38116d9e3a6SLisandro Dalcin   if (flg) {
382ce94432eSBarry Smith     if (jac->maxiter < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
383fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
38416d9e3a6SLisandro Dalcin   }
3850f1074feSSatish Balay   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_tol","Convergence tolerance PER hypre call (0.0 = use a fixed number of iterations)","None",jac->tol,&jac->tol,&flg);CHKERRQ(ierr);
38616d9e3a6SLisandro Dalcin   if (flg) {
387ce94432eSBarry Smith     if (jac->tol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be greater than or equal to zero",jac->tol);
388fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
38916d9e3a6SLisandro Dalcin   }
39016d9e3a6SLisandro Dalcin 
3910f1074feSSatish Balay   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr);
39216d9e3a6SLisandro Dalcin   if (flg) {
393ce94432eSBarry Smith     if (jac->truncfactor < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor);
394fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
39516d9e3a6SLisandro Dalcin   }
39616d9e3a6SLisandro Dalcin 
3970f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_P_max","Max elements per row for interpolation operator (0=unlimited)","None",jac->pmax,&jac->pmax,&flg);CHKERRQ(ierr);
3980f1074feSSatish Balay   if (flg) {
399ce94432eSBarry Smith     if (jac->pmax < 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"P_max %G must be greater than or equal to zero",jac->pmax);
400fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
4010f1074feSSatish Balay   }
4020f1074feSSatish Balay 
4030f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);CHKERRQ(ierr);
4040f1074feSSatish Balay   if (flg) {
405ce94432eSBarry Smith     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",jac->agg_nl);
4060f1074feSSatish Balay 
407fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
4080f1074feSSatish Balay   }
4090f1074feSSatish Balay 
4100f1074feSSatish Balay 
4110f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_num_paths","Number of paths for aggressive coarsening","None",jac->agg_num_paths,&jac->agg_num_paths,&flg);CHKERRQ(ierr);
4120f1074feSSatish Balay   if (flg) {
413ce94432eSBarry Smith     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",jac->agg_num_paths);
4140f1074feSSatish Balay 
415fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
4160f1074feSSatish Balay   }
4170f1074feSSatish Balay 
4180f1074feSSatish Balay 
41916d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr);
42016d9e3a6SLisandro Dalcin   if (flg) {
421ce94432eSBarry Smith     if (jac->strongthreshold < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold);
422fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
42316d9e3a6SLisandro Dalcin   }
42416d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr);
42516d9e3a6SLisandro Dalcin   if (flg) {
426ce94432eSBarry Smith     if (jac->maxrowsum < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum);
427ce94432eSBarry Smith     if (jac->maxrowsum > 1.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum);
428fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
42916d9e3a6SLisandro Dalcin   }
43016d9e3a6SLisandro Dalcin 
43116d9e3a6SLisandro Dalcin   /* Grid sweeps */
4320f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_all","Number of sweeps for the up and down grid levels","None",jac->gridsweeps[0],&indx,&flg);CHKERRQ(ierr);
43316d9e3a6SLisandro Dalcin   if (flg) {
434fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver,indx));
43516d9e3a6SLisandro Dalcin     /* modify the jac structure so we can view the updated options with PC_View */
43616d9e3a6SLisandro Dalcin     jac->gridsweeps[0] = indx;
4370f1074feSSatish Balay     jac->gridsweeps[1] = indx;
4380f1074feSSatish Balay     /*defaults coarse to 1 */
4390f1074feSSatish Balay     jac->gridsweeps[2] = 1;
44016d9e3a6SLisandro Dalcin   }
4410f1074feSSatish Balay 
4420f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx,&flg);CHKERRQ(ierr);
44316d9e3a6SLisandro Dalcin   if (flg) {
444fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 1));
4450f1074feSSatish Balay     jac->gridsweeps[0] = indx;
44616d9e3a6SLisandro Dalcin   }
44716d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr);
44816d9e3a6SLisandro Dalcin   if (flg) {
449fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 2));
4500f1074feSSatish Balay     jac->gridsweeps[1] = indx;
45116d9e3a6SLisandro Dalcin   }
4520f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr);
45316d9e3a6SLisandro Dalcin   if (flg) {
454fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 3));
4550f1074feSSatish Balay     jac->gridsweeps[2] = indx;
45616d9e3a6SLisandro Dalcin   }
45716d9e3a6SLisandro Dalcin 
45816d9e3a6SLisandro Dalcin   /* Relax type */
459a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_all","Relax type for the up and down cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
46016d9e3a6SLisandro Dalcin   if (flg) {
4610f1074feSSatish Balay     jac->relaxtype[0] = jac->relaxtype[1]  = indx;
462fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, indx));
4630f1074feSSatish Balay     /* by default, coarse type set to 9 */
4640f1074feSSatish Balay     jac->relaxtype[2] = 9;
4650f1074feSSatish Balay 
46616d9e3a6SLisandro Dalcin   }
467a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
46816d9e3a6SLisandro Dalcin   if (flg) {
46916d9e3a6SLisandro Dalcin     jac->relaxtype[0] = indx;
470fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 1));
47116d9e3a6SLisandro Dalcin   }
472a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
47316d9e3a6SLisandro Dalcin   if (flg) {
4740f1074feSSatish Balay     jac->relaxtype[1] = indx;
475fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 2));
47616d9e3a6SLisandro Dalcin   }
477a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr);
47816d9e3a6SLisandro Dalcin   if (flg) {
4790f1074feSSatish Balay     jac->relaxtype[2] = indx;
480fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 3));
48116d9e3a6SLisandro Dalcin   }
48216d9e3a6SLisandro Dalcin 
48316d9e3a6SLisandro Dalcin   /* Relaxation Weight */
48416d9e3a6SLisandro Dalcin   ierr = 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);CHKERRQ(ierr);
48516d9e3a6SLisandro Dalcin   if (flg) {
486fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxWt,(jac->hsolver,tmpdbl));
48716d9e3a6SLisandro Dalcin     jac->relaxweight = tmpdbl;
48816d9e3a6SLisandro Dalcin   }
48916d9e3a6SLisandro Dalcin 
49016d9e3a6SLisandro Dalcin   n         = 2;
49116d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
49216d9e3a6SLisandro Dalcin   ierr      = PetscOptionsRealArray("-pc_hypre_boomeramg_relax_weight_level","Set the relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);CHKERRQ(ierr);
49316d9e3a6SLisandro Dalcin   if (flg) {
49416d9e3a6SLisandro Dalcin     if (n == 2) {
49516d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
496fd3f9acdSBarry Smith       PetscStackCallStandard(HYPRE_BoomerAMGSetLevelRelaxWt,(jac->hsolver,twodbl[0],indx));
497ce94432eSBarry Smith     } 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);
49816d9e3a6SLisandro Dalcin   }
49916d9e3a6SLisandro Dalcin 
50016d9e3a6SLisandro Dalcin   /* Outer relaxation Weight */
50116d9e3a6SLisandro Dalcin   ierr = 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);CHKERRQ(ierr);
50216d9e3a6SLisandro Dalcin   if (flg) {
503fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetOuterWt,(jac->hsolver, tmpdbl));
50416d9e3a6SLisandro Dalcin     jac->outerrelaxweight = tmpdbl;
50516d9e3a6SLisandro Dalcin   }
50616d9e3a6SLisandro Dalcin 
50716d9e3a6SLisandro Dalcin   n         = 2;
50816d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
50916d9e3a6SLisandro Dalcin   ierr      = PetscOptionsRealArray("-pc_hypre_boomeramg_outer_relax_weight_level","Set the outer relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);CHKERRQ(ierr);
51016d9e3a6SLisandro Dalcin   if (flg) {
51116d9e3a6SLisandro Dalcin     if (n == 2) {
51216d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
513fd3f9acdSBarry Smith       PetscStackCallStandard(HYPRE_BoomerAMGSetLevelOuterWt,(jac->hsolver, twodbl[0], indx));
514ce94432eSBarry Smith     } 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);
51516d9e3a6SLisandro Dalcin   }
51616d9e3a6SLisandro Dalcin 
51716d9e3a6SLisandro Dalcin   /* the Relax Order */
518acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
51916d9e3a6SLisandro Dalcin 
52016d9e3a6SLisandro Dalcin   if (flg) {
52116d9e3a6SLisandro Dalcin     jac->relaxorder = 0;
522fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
52316d9e3a6SLisandro Dalcin   }
524a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,ALEN(HYPREBoomerAMGMeasureType),HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr);
52516d9e3a6SLisandro Dalcin   if (flg) {
52616d9e3a6SLisandro Dalcin     jac->measuretype = indx;
527fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
52816d9e3a6SLisandro Dalcin   }
5290f1074feSSatish Balay   /* update list length 3/07 */
530a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,ALEN(HYPREBoomerAMGCoarsenType),HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr);
53116d9e3a6SLisandro Dalcin   if (flg) {
53216d9e3a6SLisandro Dalcin     jac->coarsentype = indx;
533fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
53416d9e3a6SLisandro Dalcin   }
5350f1074feSSatish Balay 
5360f1074feSSatish Balay   /* new 3/07 */
537a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,ALEN(HYPREBoomerAMGInterpType),HYPREBoomerAMGInterpType[0],&indx,&flg);CHKERRQ(ierr);
5380f1074feSSatish Balay   if (flg) {
5390f1074feSSatish Balay     jac->interptype = indx;
540fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
5410f1074feSSatish Balay   }
5420f1074feSSatish Balay 
543b96a4a96SBarry Smith   ierr = PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);CHKERRQ(ierr);
54416d9e3a6SLisandro Dalcin   if (flg) {
545b96a4a96SBarry Smith     level = 3;
5460298fd71SBarry Smith     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,NULL);CHKERRQ(ierr);
5472fa5cd67SKarl Rupp 
548b96a4a96SBarry Smith     jac->printstatistics = PETSC_TRUE;
549fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetPrintLevel,(jac->hsolver,level));
5502ae77aedSBarry Smith   }
5512ae77aedSBarry Smith 
552b96a4a96SBarry Smith   ierr = PetscOptionsName("-pc_hypre_boomeramg_print_debug","Print debug information","None",&flg);CHKERRQ(ierr);
5532ae77aedSBarry Smith   if (flg) {
554b96a4a96SBarry Smith     level = 3;
5550298fd71SBarry Smith     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,NULL);CHKERRQ(ierr);
5562fa5cd67SKarl Rupp 
557b96a4a96SBarry Smith     jac->printstatistics = PETSC_TRUE;
558fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetDebugFlag,(jac->hsolver,level));
55916d9e3a6SLisandro Dalcin   }
5608f87f92bSBarry Smith 
561acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5628f87f92bSBarry Smith   if (flg && tmp_truth) {
5638f87f92bSBarry Smith     jac->nodal_coarsen = 1;
564fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetNodal,(jac->hsolver,1));
5658f87f92bSBarry Smith   }
5668f87f92bSBarry Smith 
567acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5688f87f92bSBarry Smith   if (flg && tmp_truth) {
5698f87f92bSBarry Smith     PetscInt tmp_int;
5708f87f92bSBarry Smith     ierr = PetscOptionsInt("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);CHKERRQ(ierr);
5718f87f92bSBarry Smith     if (flg) jac->nodal_relax_levels = tmp_int;
572fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothType,(jac->hsolver,6));
573fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetDomainType,(jac->hsolver,1));
574fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetOverlap,(jac->hsolver,0));
575fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothNumLevels,(jac->hsolver,jac->nodal_relax_levels));
5768f87f92bSBarry Smith   }
5778f87f92bSBarry Smith 
57816d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
57916d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
58016d9e3a6SLisandro Dalcin }
58116d9e3a6SLisandro Dalcin 
58216d9e3a6SLisandro Dalcin #undef __FUNCT__
58316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG"
584ace3abfcSBarry Smith 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)
58516d9e3a6SLisandro Dalcin {
58616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
58716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
5884ddd07fcSJed Brown   PetscInt       oits;
58916d9e3a6SLisandro Dalcin 
59016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
591fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,its*jac->maxiter));
592fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,rtol));
59316d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_TRUE;
59416d9e3a6SLisandro Dalcin   ierr                 = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr);
59516d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_FALSE;
596fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGGetNumIterations,(jac->hsolver,&oits));
5974d0a8057SBarry Smith   *outits = oits;
5984d0a8057SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
5994d0a8057SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
600fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
601fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
60216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
60316d9e3a6SLisandro Dalcin }
60416d9e3a6SLisandro Dalcin 
60516d9e3a6SLisandro Dalcin 
60616d9e3a6SLisandro Dalcin #undef __FUNCT__
60716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_BoomerAMG"
60816d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
60916d9e3a6SLisandro Dalcin {
61016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
61116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
612ace3abfcSBarry Smith   PetscBool      iascii;
61316d9e3a6SLisandro Dalcin 
61416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
615251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
61616d9e3a6SLisandro Dalcin   if (iascii) {
61716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr);
61816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr);
61916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr);
62016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr);
62116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr);
62216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr);
6230f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation truncation factor %G\n",jac->truncfactor);CHKERRQ(ierr);
6240f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);CHKERRQ(ierr);
6250f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);CHKERRQ(ierr);
6260f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);CHKERRQ(ierr);
6270f1074feSSatish Balay 
62816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr);
62916d9e3a6SLisandro Dalcin 
6300f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps down         %d\n",jac->gridsweeps[0]);CHKERRQ(ierr);
6310f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps up           %d\n",jac->gridsweeps[1]);CHKERRQ(ierr);
6320f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on coarse    %d\n",jac->gridsweeps[2]);CHKERRQ(ierr);
63316d9e3a6SLisandro Dalcin 
6340f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr);
6350f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr);
6360f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr);
63716d9e3a6SLisandro Dalcin 
63816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax weight  (all)      %G\n",jac->relaxweight);CHKERRQ(ierr);
63916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr);
64016d9e3a6SLisandro Dalcin 
64116d9e3a6SLisandro Dalcin     if (jac->relaxorder) {
64216d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr);
64316d9e3a6SLisandro Dalcin     } else {
64416d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr);
64516d9e3a6SLisandro Dalcin     }
64616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr);
64716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr);
6480f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation type  %s\n",HYPREBoomerAMGInterpType[jac->interptype]);CHKERRQ(ierr);
6498f87f92bSBarry Smith     if (jac->nodal_coarsen) {
6508f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal())\n");CHKERRQ(ierr);
6518f87f92bSBarry Smith     }
6528f87f92bSBarry Smith     if (jac->nodal_relax) {
6538f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);CHKERRQ(ierr);
6548f87f92bSBarry Smith     }
65516d9e3a6SLisandro Dalcin   }
65616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
65716d9e3a6SLisandro Dalcin }
65816d9e3a6SLisandro Dalcin 
65916d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
66016d9e3a6SLisandro Dalcin #undef __FUNCT__
66116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails"
66216d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc)
66316d9e3a6SLisandro Dalcin {
66416d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
66516d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
6664ddd07fcSJed Brown   PetscInt       indx;
667ace3abfcSBarry Smith   PetscBool      flag;
66816d9e3a6SLisandro Dalcin   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};
66916d9e3a6SLisandro Dalcin 
67016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
67116d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr);
67216d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr);
67316d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr);
6742fa5cd67SKarl Rupp   if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels));
67516d9e3a6SLisandro Dalcin 
67616d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr);
6772fa5cd67SKarl Rupp   if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter));
67816d9e3a6SLisandro Dalcin 
67916d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr);
6802fa5cd67SKarl Rupp   if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal));
68116d9e3a6SLisandro Dalcin 
682acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_parasails_logging","Print info to screen","None",(PetscBool)jac->logging,(PetscBool*)&jac->logging,&flag);CHKERRQ(ierr);
6832fa5cd67SKarl Rupp   if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging));
68416d9e3a6SLisandro Dalcin 
685acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscBool)jac->ruse,(PetscBool*)&jac->ruse,&flag);CHKERRQ(ierr);
6862fa5cd67SKarl Rupp   if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse));
68716d9e3a6SLisandro Dalcin 
688a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,ALEN(symtlist),symtlist[0],&indx,&flag);CHKERRQ(ierr);
68916d9e3a6SLisandro Dalcin   if (flag) {
69016d9e3a6SLisandro Dalcin     jac->symt = indx;
691fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
69216d9e3a6SLisandro Dalcin   }
69316d9e3a6SLisandro Dalcin 
69416d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
69516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
69616d9e3a6SLisandro Dalcin }
69716d9e3a6SLisandro Dalcin 
69816d9e3a6SLisandro Dalcin #undef __FUNCT__
69916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_ParaSails"
70016d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
70116d9e3a6SLisandro Dalcin {
70216d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
70316d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
704ace3abfcSBarry Smith   PetscBool      iascii;
70516d9e3a6SLisandro Dalcin   const char     *symt = 0;;
70616d9e3a6SLisandro Dalcin 
70716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
708251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
70916d9e3a6SLisandro Dalcin   if (iascii) {
71016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");CHKERRQ(ierr);
71116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr);
71216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr);
71316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr);
71416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr);
715ace3abfcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscBools[jac->ruse]);CHKERRQ(ierr);
716ace3abfcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscBools[jac->logging]);CHKERRQ(ierr);
7172fa5cd67SKarl Rupp     if (!jac->symt) symt = "nonsymmetric matrix and preconditioner";
7182fa5cd67SKarl Rupp     else if (jac->symt == 1) symt = "SPD matrix and preconditioner";
7192fa5cd67SKarl Rupp     else if (jac->symt == 2) symt = "nonsymmetric matrix but SPD preconditioner";
720ce94432eSBarry Smith     else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
72116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr);
72216d9e3a6SLisandro Dalcin   }
72316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
72416d9e3a6SLisandro Dalcin }
72516d9e3a6SLisandro Dalcin /* ---------------------------------------------------------------------------------*/
72616d9e3a6SLisandro Dalcin 
72716d9e3a6SLisandro Dalcin #undef __FUNCT__
72816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType_HYPRE"
729f7a08781SBarry Smith static PetscErrorCode  PCHYPREGetType_HYPRE(PC pc,const char *name[])
73016d9e3a6SLisandro Dalcin {
73116d9e3a6SLisandro Dalcin   PC_HYPRE *jac = (PC_HYPRE*)pc->data;
73216d9e3a6SLisandro Dalcin 
73316d9e3a6SLisandro Dalcin   PetscFunctionBegin;
73416d9e3a6SLisandro Dalcin   *name = jac->hypre_type;
73516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
73616d9e3a6SLisandro Dalcin }
73716d9e3a6SLisandro Dalcin 
73816d9e3a6SLisandro Dalcin #undef __FUNCT__
73916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType_HYPRE"
740f7a08781SBarry Smith static PetscErrorCode  PCHYPRESetType_HYPRE(PC pc,const char name[])
74116d9e3a6SLisandro Dalcin {
74216d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
74316d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
744ace3abfcSBarry Smith   PetscBool      flag;
74516d9e3a6SLisandro Dalcin 
74616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
74716d9e3a6SLisandro Dalcin   if (jac->hypre_type) {
74816d9e3a6SLisandro Dalcin     ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr);
749ce94432eSBarry Smith     if (!flag) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
75016d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
75116d9e3a6SLisandro Dalcin   } else {
75216d9e3a6SLisandro Dalcin     ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr);
75316d9e3a6SLisandro Dalcin   }
75416d9e3a6SLisandro Dalcin 
75516d9e3a6SLisandro Dalcin   jac->maxiter         = PETSC_DEFAULT;
75616d9e3a6SLisandro Dalcin   jac->tol             = PETSC_DEFAULT;
75716d9e3a6SLisandro Dalcin   jac->printstatistics = PetscLogPrintInfo;
75816d9e3a6SLisandro Dalcin 
75916d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr);
76016d9e3a6SLisandro Dalcin   if (flag) {
761fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParCSRPilutCreate,(jac->comm_hypre,&jac->hsolver));
76216d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
76316d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Pilut;
76416d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParCSRPilutDestroy;
76516d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParCSRPilutSetup;
76616d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParCSRPilutSolve;
76716d9e3a6SLisandro Dalcin     jac->factorrowsize      = PETSC_DEFAULT;
76816d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
76916d9e3a6SLisandro Dalcin   }
77016d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr);
77116d9e3a6SLisandro Dalcin   if (flag) {
772fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsCreate,(jac->comm_hypre,&jac->hsolver));
77316d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
77416d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_ParaSails;
77516d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParaSailsDestroy;
77616d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParaSailsSetup;
77716d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParaSailsSolve;
77816d9e3a6SLisandro Dalcin     /* initialize */
77916d9e3a6SLisandro Dalcin     jac->nlevels    = 1;
78016d9e3a6SLisandro Dalcin     jac->threshhold = .1;
78116d9e3a6SLisandro Dalcin     jac->filter     = .1;
78216d9e3a6SLisandro Dalcin     jac->loadbal    = 0;
7832fa5cd67SKarl Rupp     if (PetscLogPrintInfo) jac->logging = (int) PETSC_TRUE;
7842fa5cd67SKarl Rupp     else jac->logging = (int) PETSC_FALSE;
7852fa5cd67SKarl Rupp 
78616d9e3a6SLisandro Dalcin     jac->ruse = (int) PETSC_FALSE;
78716d9e3a6SLisandro Dalcin     jac->symt = 0;
788fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels));
789fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter));
790fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal));
791fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging));
792fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse));
793fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
79416d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
79516d9e3a6SLisandro Dalcin   }
79616d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr);
79716d9e3a6SLisandro Dalcin   if (flag) {
79816d9e3a6SLisandro Dalcin     ierr                    = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver);
79916d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
80016d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Euclid;
80116d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_EuclidDestroy;
80216d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_EuclidSetup;
80316d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_EuclidSolve;
80416d9e3a6SLisandro Dalcin     /* initialization */
80516d9e3a6SLisandro Dalcin     jac->bjilu              = PETSC_FALSE;
80616d9e3a6SLisandro Dalcin     jac->levels             = 1;
80716d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
80816d9e3a6SLisandro Dalcin   }
80916d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr);
81016d9e3a6SLisandro Dalcin   if (flag) {
81116d9e3a6SLisandro Dalcin     ierr                     = HYPRE_BoomerAMGCreate(&jac->hsolver);
81216d9e3a6SLisandro Dalcin     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
81316d9e3a6SLisandro Dalcin     pc->ops->view            = PCView_HYPRE_BoomerAMG;
81416d9e3a6SLisandro Dalcin     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
81516d9e3a6SLisandro Dalcin     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
81616d9e3a6SLisandro Dalcin     jac->destroy             = HYPRE_BoomerAMGDestroy;
81716d9e3a6SLisandro Dalcin     jac->setup               = HYPRE_BoomerAMGSetup;
81816d9e3a6SLisandro Dalcin     jac->solve               = HYPRE_BoomerAMGSolve;
81916d9e3a6SLisandro Dalcin     jac->applyrichardson     = PETSC_FALSE;
82016d9e3a6SLisandro Dalcin     /* these defaults match the hypre defaults */
82116d9e3a6SLisandro Dalcin     jac->cycletype        = 1;
82216d9e3a6SLisandro Dalcin     jac->maxlevels        = 25;
82316d9e3a6SLisandro Dalcin     jac->maxiter          = 1;
8248f87f92bSBarry Smith     jac->tol              = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */
82516d9e3a6SLisandro Dalcin     jac->truncfactor      = 0.0;
82616d9e3a6SLisandro Dalcin     jac->strongthreshold  = .25;
82716d9e3a6SLisandro Dalcin     jac->maxrowsum        = .9;
82816d9e3a6SLisandro Dalcin     jac->coarsentype      = 6;
82916d9e3a6SLisandro Dalcin     jac->measuretype      = 0;
8300f1074feSSatish Balay     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = 1;
8318f87f92bSBarry Smith     jac->relaxtype[0]     = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */
8320f1074feSSatish Balay     jac->relaxtype[2]     = 9; /*G.E. */
83316d9e3a6SLisandro Dalcin     jac->relaxweight      = 1.0;
83416d9e3a6SLisandro Dalcin     jac->outerrelaxweight = 1.0;
83516d9e3a6SLisandro Dalcin     jac->relaxorder       = 1;
8360f1074feSSatish Balay     jac->interptype       = 0;
8370f1074feSSatish Balay     jac->agg_nl           = 0;
8380f1074feSSatish Balay     jac->pmax             = 0;
8390f1074feSSatish Balay     jac->truncfactor      = 0.0;
8400f1074feSSatish Balay     jac->agg_num_paths    = 1;
8418f87f92bSBarry Smith 
8428f87f92bSBarry Smith     jac->nodal_coarsen      = 0;
8438f87f92bSBarry Smith     jac->nodal_relax        = PETSC_FALSE;
8448f87f92bSBarry Smith     jac->nodal_relax_levels = 1;
845fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
846fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
847fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
848fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
849fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
850fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
851fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
852fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
853fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
854fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
855fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
856fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
857fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
858fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
859fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, jac->relaxtype[0]));  /*defaults coarse to 9*/
860fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver, jac->gridsweeps[0])); /*defaults coarse to 1 */
86116d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
86216d9e3a6SLisandro Dalcin   }
863503cfb0cSBarry Smith   ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr);
8642fa5cd67SKarl Rupp 
8650298fd71SBarry Smith   jac->hypre_type = NULL;
866ce94432eSBarry Smith   SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
86716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
86816d9e3a6SLisandro Dalcin }
86916d9e3a6SLisandro Dalcin 
87016d9e3a6SLisandro Dalcin /*
87116d9e3a6SLisandro Dalcin     It only gets here if the HYPRE type has not been set before the call to
87216d9e3a6SLisandro Dalcin    ...SetFromOptions() which actually is most of the time
87316d9e3a6SLisandro Dalcin */
87416d9e3a6SLisandro Dalcin #undef __FUNCT__
87516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE"
87616d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE(PC pc)
87716d9e3a6SLisandro Dalcin {
87816d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
8794ddd07fcSJed Brown   PetscInt       indx;
88016d9e3a6SLisandro Dalcin   const char     *type[] = {"pilut","parasails","boomeramg","euclid"};
881ace3abfcSBarry Smith   PetscBool      flg;
88216d9e3a6SLisandro Dalcin 
88316d9e3a6SLisandro Dalcin   PetscFunctionBegin;
88416d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr);
88502a17cd4SBarry Smith   ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"boomeramg",&indx,&flg);CHKERRQ(ierr);
88616d9e3a6SLisandro Dalcin   if (flg) {
88716d9e3a6SLisandro Dalcin     ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr);
88802a17cd4SBarry Smith   } else {
88902a17cd4SBarry Smith     ierr = PCHYPRESetType_HYPRE(pc,"boomeramg");CHKERRQ(ierr);
89016d9e3a6SLisandro Dalcin   }
89116d9e3a6SLisandro Dalcin   if (pc->ops->setfromoptions) {
89216d9e3a6SLisandro Dalcin     ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr);
89316d9e3a6SLisandro Dalcin   }
89416d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
89516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
89616d9e3a6SLisandro Dalcin }
89716d9e3a6SLisandro Dalcin 
89816d9e3a6SLisandro Dalcin #undef __FUNCT__
89916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType"
90016d9e3a6SLisandro Dalcin /*@C
90116d9e3a6SLisandro Dalcin      PCHYPRESetType - Sets which hypre preconditioner you wish to use
90216d9e3a6SLisandro Dalcin 
90316d9e3a6SLisandro Dalcin    Input Parameters:
90416d9e3a6SLisandro Dalcin +     pc - the preconditioner context
90516d9e3a6SLisandro Dalcin -     name - either  pilut, parasails, boomeramg, euclid
90616d9e3a6SLisandro Dalcin 
90716d9e3a6SLisandro Dalcin    Options Database Keys:
90816d9e3a6SLisandro Dalcin    -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
90916d9e3a6SLisandro Dalcin 
91016d9e3a6SLisandro Dalcin    Level: intermediate
91116d9e3a6SLisandro Dalcin 
91216d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
91316d9e3a6SLisandro Dalcin            PCHYPRE
91416d9e3a6SLisandro Dalcin 
91516d9e3a6SLisandro Dalcin @*/
9167087cfbeSBarry Smith PetscErrorCode  PCHYPRESetType(PC pc,const char name[])
91716d9e3a6SLisandro Dalcin {
9184ac538c5SBarry Smith   PetscErrorCode ierr;
91916d9e3a6SLisandro Dalcin 
92016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9210700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
92216d9e3a6SLisandro Dalcin   PetscValidCharPointer(name,2);
9234ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCHYPRESetType_C",(PC,const char[]),(pc,name));CHKERRQ(ierr);
92416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
92516d9e3a6SLisandro Dalcin }
92616d9e3a6SLisandro Dalcin 
92716d9e3a6SLisandro Dalcin #undef __FUNCT__
92816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType"
92916d9e3a6SLisandro Dalcin /*@C
93016d9e3a6SLisandro Dalcin      PCHYPREGetType - Gets which hypre preconditioner you are using
93116d9e3a6SLisandro Dalcin 
93216d9e3a6SLisandro Dalcin    Input Parameter:
93316d9e3a6SLisandro Dalcin .     pc - the preconditioner context
93416d9e3a6SLisandro Dalcin 
93516d9e3a6SLisandro Dalcin    Output Parameter:
93616d9e3a6SLisandro Dalcin .     name - either  pilut, parasails, boomeramg, euclid
93716d9e3a6SLisandro Dalcin 
93816d9e3a6SLisandro Dalcin    Level: intermediate
93916d9e3a6SLisandro Dalcin 
94016d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC,
94116d9e3a6SLisandro Dalcin            PCHYPRE
94216d9e3a6SLisandro Dalcin 
94316d9e3a6SLisandro Dalcin @*/
9447087cfbeSBarry Smith PetscErrorCode  PCHYPREGetType(PC pc,const char *name[])
94516d9e3a6SLisandro Dalcin {
9464ac538c5SBarry Smith   PetscErrorCode ierr;
94716d9e3a6SLisandro Dalcin 
94816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9490700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
95016d9e3a6SLisandro Dalcin   PetscValidPointer(name,2);
9514ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCHYPREGetType_C",(PC,const char*[]),(pc,name));CHKERRQ(ierr);
95216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
95316d9e3a6SLisandro Dalcin }
95416d9e3a6SLisandro Dalcin 
95516d9e3a6SLisandro Dalcin /*MC
95616d9e3a6SLisandro Dalcin      PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre
95716d9e3a6SLisandro Dalcin 
95816d9e3a6SLisandro Dalcin    Options Database Keys:
95916d9e3a6SLisandro Dalcin +   -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
96016d9e3a6SLisandro Dalcin -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
96116d9e3a6SLisandro Dalcin           preconditioner
96216d9e3a6SLisandro Dalcin 
96316d9e3a6SLisandro Dalcin    Level: intermediate
96416d9e3a6SLisandro Dalcin 
96516d9e3a6SLisandro Dalcin    Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()),
96616d9e3a6SLisandro Dalcin           the many hypre options can ONLY be set via the options database (e.g. the command line
96716d9e3a6SLisandro Dalcin           or with PetscOptionsSetValue(), there are no functions to set them)
96816d9e3a6SLisandro Dalcin 
96916d9e3a6SLisandro Dalcin           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations
9700f1074feSSatish Balay           (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if
9710f1074feSSatish Balay           -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner
9720f1074feSSatish Balay           (-pc_hypre_boomeramg_rtol should be set to 0.0 - the default - to strictly use a fixed number of
9738f87f92bSBarry Smith           iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations
9740f1074feSSatish Balay           and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10
9750f1074feSSatish Balay           then AT MOST twenty V-cycles of boomeramg will be called.
97616d9e3a6SLisandro Dalcin 
9770f1074feSSatish Balay            Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation
9780f1074feSSatish Balay            (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry.
9790f1074feSSatish Balay            Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi.
98016d9e3a6SLisandro Dalcin           If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly
98116d9e3a6SLisandro Dalcin           and use -ksp_max_it to control the number of V-cycles.
98216d9e3a6SLisandro Dalcin           (see the PETSc FAQ.html at the PETSc website under the Documentation tab).
98316d9e3a6SLisandro Dalcin 
98416d9e3a6SLisandro Dalcin           2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option
98516d9e3a6SLisandro Dalcin           -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L.
98616d9e3a6SLisandro Dalcin 
9879e5bc791SBarry Smith           See PCPFMG for access to the hypre Struct PFMG solver
9889e5bc791SBarry Smith 
98916d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
9909e5bc791SBarry Smith            PCHYPRESetType(), PCPFMG
99116d9e3a6SLisandro Dalcin 
99216d9e3a6SLisandro Dalcin M*/
99316d9e3a6SLisandro Dalcin 
99416d9e3a6SLisandro Dalcin #undef __FUNCT__
99516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCCreate_HYPRE"
996*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_HYPRE(PC pc)
99716d9e3a6SLisandro Dalcin {
99816d9e3a6SLisandro Dalcin   PC_HYPRE       *jac;
99916d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
100016d9e3a6SLisandro Dalcin 
100116d9e3a6SLisandro Dalcin   PetscFunctionBegin;
100238f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_HYPRE,&jac);CHKERRQ(ierr);
10032fa5cd67SKarl Rupp 
100416d9e3a6SLisandro Dalcin   pc->data                = jac;
100516d9e3a6SLisandro Dalcin   pc->ops->destroy        = PCDestroy_HYPRE;
100616d9e3a6SLisandro Dalcin   pc->ops->setfromoptions = PCSetFromOptions_HYPRE;
100716d9e3a6SLisandro Dalcin   pc->ops->setup          = PCSetUp_HYPRE;
100816d9e3a6SLisandro Dalcin   pc->ops->apply          = PCApply_HYPRE;
100916d9e3a6SLisandro Dalcin   jac->comm_hypre         = MPI_COMM_NULL;
10100298fd71SBarry Smith   jac->hypre_type         = NULL;
101116d9e3a6SLisandro Dalcin   /* duplicate communicator for hypre */
1012ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(jac->comm_hypre));CHKERRQ(ierr);
101300de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr);
101400de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr);
101516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
101616d9e3a6SLisandro Dalcin }
1017ebc551c0SBarry Smith 
1018f91d8e95SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------*/
1019f91d8e95SBarry Smith 
1020b862ddfaSBarry Smith /* this include is needed ONLY to allow access to the private data inside the Mat object specific to hypre */
1021b45d2f2cSJed Brown #include <petsc-private/matimpl.h>
1022ebc551c0SBarry Smith 
1023ebc551c0SBarry Smith typedef struct {
102468326731SBarry Smith   MPI_Comm           hcomm;        /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */
1025f91d8e95SBarry Smith   HYPRE_StructSolver hsolver;
10269e5bc791SBarry Smith 
10279e5bc791SBarry Smith   /* keep copy of PFMG options used so may view them */
10284ddd07fcSJed Brown   PetscInt its;
10299e5bc791SBarry Smith   double   tol;
10304ddd07fcSJed Brown   PetscInt relax_type;
10314ddd07fcSJed Brown   PetscInt rap_type;
10324ddd07fcSJed Brown   PetscInt num_pre_relax,num_post_relax;
10334ddd07fcSJed Brown   PetscInt max_levels;
1034ebc551c0SBarry Smith } PC_PFMG;
1035ebc551c0SBarry Smith 
1036ebc551c0SBarry Smith #undef __FUNCT__
1037ebc551c0SBarry Smith #define __FUNCT__ "PCDestroy_PFMG"
1038ebc551c0SBarry Smith PetscErrorCode PCDestroy_PFMG(PC pc)
1039ebc551c0SBarry Smith {
1040ebc551c0SBarry Smith   PetscErrorCode ierr;
1041f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1042ebc551c0SBarry Smith 
1043ebc551c0SBarry Smith   PetscFunctionBegin;
10442fa5cd67SKarl Rupp   if (ex->hsolver) PetscStackCallStandard(HYPRE_StructPFMGDestroy,(ex->hsolver));
1045f91d8e95SBarry Smith   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1046c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1047ebc551c0SBarry Smith   PetscFunctionReturn(0);
1048ebc551c0SBarry Smith }
1049ebc551c0SBarry Smith 
10509e5bc791SBarry Smith static const char *PFMGRelaxType[] = {"Jacobi","Weighted-Jacobi","symmetric-Red/Black-Gauss-Seidel","Red/Black-Gauss-Seidel"};
10519e5bc791SBarry Smith static const char *PFMGRAPType[] = {"Galerkin","non-Galerkin"};
10529e5bc791SBarry Smith 
1053ebc551c0SBarry Smith #undef __FUNCT__
1054ebc551c0SBarry Smith #define __FUNCT__ "PCView_PFMG"
1055ebc551c0SBarry Smith PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer)
1056ebc551c0SBarry Smith {
1057ebc551c0SBarry Smith   PetscErrorCode ierr;
1058ace3abfcSBarry Smith   PetscBool      iascii;
1059f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1060ebc551c0SBarry Smith 
1061ebc551c0SBarry Smith   PetscFunctionBegin;
1062251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10639e5bc791SBarry Smith   if (iascii) {
10649e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG preconditioning\n");CHKERRQ(ierr);
10659e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
10669e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
10679e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
10689e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: RAP type %s\n",PFMGRAPType[ex->rap_type]);CHKERRQ(ierr);
10699e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
10703b46a515SGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max levels %d\n",ex->max_levels);CHKERRQ(ierr);
10719e5bc791SBarry Smith   }
1072ebc551c0SBarry Smith   PetscFunctionReturn(0);
1073ebc551c0SBarry Smith }
1074ebc551c0SBarry Smith 
10759e5bc791SBarry Smith 
1076ebc551c0SBarry Smith #undef __FUNCT__
1077ebc551c0SBarry Smith #define __FUNCT__ "PCSetFromOptions_PFMG"
1078ebc551c0SBarry Smith PetscErrorCode PCSetFromOptions_PFMG(PC pc)
1079ebc551c0SBarry Smith {
1080ebc551c0SBarry Smith   PetscErrorCode ierr;
1081f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1082ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
1083ebc551c0SBarry Smith 
1084ebc551c0SBarry Smith   PetscFunctionBegin;
1085ebc551c0SBarry Smith   ierr = PetscOptionsHead("PFMG options");CHKERRQ(ierr);
10860298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",flg,&flg,NULL);CHKERRQ(ierr);
108768326731SBarry Smith   if (flg) {
1088a0324ebeSBarry Smith     int level=3;
1089fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_StructPFMGSetPrintLevel,(ex->hsolver,level));
109068326731SBarry Smith   }
10910298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_pfmg_its","Number of iterations of PFMG to use as preconditioner","HYPRE_StructPFMGSetMaxIter",ex->its,&ex->its,NULL);CHKERRQ(ierr);
1092fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(ex->hsolver,ex->its));
10930298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_pfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_StructPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,NULL);CHKERRQ(ierr);
1094fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetNumPreRelax,(ex->hsolver,ex->num_pre_relax));
10950298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_pfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_StructPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,NULL);CHKERRQ(ierr);
1096fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetNumPostRelax,(ex->hsolver,ex->num_post_relax));
10979e5bc791SBarry Smith 
10980298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_pfmg_max_levels","Max Levels for MG hierarchy","HYPRE_StructPFMGSetMaxLevels",ex->max_levels,&ex->max_levels,NULL);CHKERRQ(ierr);
1099fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetMaxLevels,(ex->hsolver,ex->max_levels));
11003b46a515SGlenn Hammond 
11010298fd71SBarry Smith   ierr = PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,NULL);CHKERRQ(ierr);
1102fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(ex->hsolver,ex->tol));
11030298fd71SBarry Smith   ierr = 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);CHKERRQ(ierr);
1104fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetRelaxType,(ex->hsolver, ex->relax_type));
11050298fd71SBarry Smith   ierr = PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,ALEN(PFMGRAPType),PFMGRAPType[ex->rap_type],&ex->rap_type,NULL);CHKERRQ(ierr);
1106fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetRAPType,(ex->hsolver, ex->rap_type));
1107ebc551c0SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1108ebc551c0SBarry Smith   PetscFunctionReturn(0);
1109ebc551c0SBarry Smith }
1110ebc551c0SBarry Smith 
1111f91d8e95SBarry Smith #undef __FUNCT__
1112f91d8e95SBarry Smith #define __FUNCT__ "PCApply_PFMG"
1113f91d8e95SBarry Smith PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y)
1114f91d8e95SBarry Smith {
1115f91d8e95SBarry Smith   PetscErrorCode  ierr;
1116f91d8e95SBarry Smith   PC_PFMG         *ex = (PC_PFMG*) pc->data;
1117f91d8e95SBarry Smith   PetscScalar     *xx,*yy;
11184ddd07fcSJed Brown   PetscInt        ilower[3],iupper[3];
111968326731SBarry Smith   Mat_HYPREStruct *mx = (Mat_HYPREStruct*)(pc->pmat->data);
1120f91d8e95SBarry Smith 
1121f91d8e95SBarry Smith   PetscFunctionBegin;
1122aa219208SBarry Smith   ierr       = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1123f91d8e95SBarry Smith   iupper[0] += ilower[0] - 1;
1124f91d8e95SBarry Smith   iupper[1] += ilower[1] - 1;
1125f91d8e95SBarry Smith   iupper[2] += ilower[2] - 1;
1126f91d8e95SBarry Smith 
1127f91d8e95SBarry Smith   /* copy x values over to hypre */
1128fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorSetConstantValues,(mx->hb,0.0));
1129f91d8e95SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1130fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorSetBoxValues,(mx->hb,ilower,iupper,xx));
1131f91d8e95SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1132fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorAssemble,(mx->hb));
1133fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSolve,(ex->hsolver,mx->hmat,mx->hb,mx->hx));
1134f91d8e95SBarry Smith 
1135f91d8e95SBarry Smith   /* copy solution values back to PETSc */
1136f91d8e95SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1137fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorGetBoxValues,(mx->hx,ilower,iupper,yy));
1138f91d8e95SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1139f91d8e95SBarry Smith   PetscFunctionReturn(0);
1140f91d8e95SBarry Smith }
1141f91d8e95SBarry Smith 
11429e5bc791SBarry Smith #undef __FUNCT__
11439e5bc791SBarry Smith #define __FUNCT__ "PCApplyRichardson_PFMG"
1144ace3abfcSBarry Smith 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)
11459e5bc791SBarry Smith {
11469e5bc791SBarry Smith   PC_PFMG        *jac = (PC_PFMG*)pc->data;
11479e5bc791SBarry Smith   PetscErrorCode ierr;
11484ddd07fcSJed Brown   PetscInt       oits;
11499e5bc791SBarry Smith 
11509e5bc791SBarry Smith   PetscFunctionBegin;
1151fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(jac->hsolver,its*jac->its));
1152fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(jac->hsolver,rtol));
11539e5bc791SBarry Smith 
11549e5bc791SBarry Smith   ierr = PCApply_PFMG(pc,b,y);CHKERRQ(ierr);
1155fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGGetNumIterations,(jac->hsolver,&oits));
11569e5bc791SBarry Smith   *outits = oits;
11579e5bc791SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
11589e5bc791SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
1159fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(jac->hsolver,jac->tol));
1160fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(jac->hsolver,jac->its));
11619e5bc791SBarry Smith   PetscFunctionReturn(0);
11629e5bc791SBarry Smith }
11639e5bc791SBarry Smith 
11649e5bc791SBarry Smith 
11653a32d3dbSGlenn Hammond #undef __FUNCT__
11663a32d3dbSGlenn Hammond #define __FUNCT__ "PCSetUp_PFMG"
11673a32d3dbSGlenn Hammond PetscErrorCode PCSetUp_PFMG(PC pc)
11683a32d3dbSGlenn Hammond {
11693a32d3dbSGlenn Hammond   PetscErrorCode  ierr;
11703a32d3dbSGlenn Hammond   PC_PFMG         *ex = (PC_PFMG*) pc->data;
11713a32d3dbSGlenn Hammond   Mat_HYPREStruct *mx = (Mat_HYPREStruct*)(pc->pmat->data);
1172ace3abfcSBarry Smith   PetscBool       flg;
11733a32d3dbSGlenn Hammond 
11743a32d3dbSGlenn Hammond   PetscFunctionBegin;
1175251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESTRUCT,&flg);CHKERRQ(ierr);
1176ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESTRUCT with this preconditioner");
11773a32d3dbSGlenn Hammond 
11783a32d3dbSGlenn Hammond   /* create the hypre solver object and set its information */
11792fa5cd67SKarl Rupp   if (ex->hsolver) PetscStackCallStandard(HYPRE_StructPFMGDestroy,(ex->hsolver));
1180fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
11813a32d3dbSGlenn Hammond   ierr = PCSetFromOptions_PFMG(pc);CHKERRQ(ierr);
1182fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetup,(ex->hsolver,mx->hmat,mx->hb,mx->hx));
1183fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetZeroGuess,(ex->hsolver));
11843a32d3dbSGlenn Hammond   PetscFunctionReturn(0);
11853a32d3dbSGlenn Hammond }
11863a32d3dbSGlenn Hammond 
1187ebc551c0SBarry Smith 
1188ebc551c0SBarry Smith /*MC
1189ebc551c0SBarry Smith      PCPFMG - the hypre PFMG multigrid solver
1190ebc551c0SBarry Smith 
1191ebc551c0SBarry Smith    Level: advanced
1192ebc551c0SBarry Smith 
11939e5bc791SBarry Smith    Options Database:
11949e5bc791SBarry Smith + -pc_pfmg_its <its> number of iterations of PFMG to use as preconditioner
11959e5bc791SBarry Smith . -pc_pfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
11969e5bc791SBarry Smith . -pc_pfmg_num_post_relax <steps> number of smoothing steps after coarse grid
11979e5bc791SBarry Smith . -pc_pfmg_tol <tol> tolerance of PFMG
11989e5bc791SBarry Smith . -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
11999e5bc791SBarry Smith - -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin
1200f91d8e95SBarry Smith 
12019e5bc791SBarry Smith    Notes:  This is for CELL-centered descretizations
12029e5bc791SBarry Smith 
12038e395302SJed Brown            This must be used with the MATHYPRESTRUCT matrix type.
1204aa219208SBarry Smith            This is less general than in hypre, it supports only one block per process defined by a PETSc DMDA.
12059e5bc791SBarry Smith 
12069e5bc791SBarry Smith .seealso:  PCMG, MATHYPRESTRUCT
1207ebc551c0SBarry Smith M*/
1208ebc551c0SBarry Smith 
1209ebc551c0SBarry Smith #undef __FUNCT__
1210ebc551c0SBarry Smith #define __FUNCT__ "PCCreate_PFMG"
1211*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_PFMG(PC pc)
1212ebc551c0SBarry Smith {
1213ebc551c0SBarry Smith   PetscErrorCode ierr;
1214ebc551c0SBarry Smith   PC_PFMG        *ex;
1215ebc551c0SBarry Smith 
1216ebc551c0SBarry Smith   PetscFunctionBegin;
1217ebc551c0SBarry Smith   ierr     = PetscNew(PC_PFMG,&ex);CHKERRQ(ierr); \
121868326731SBarry Smith   pc->data = ex;
1219ebc551c0SBarry Smith 
12209e5bc791SBarry Smith   ex->its            = 1;
12219e5bc791SBarry Smith   ex->tol            = 1.e-8;
12229e5bc791SBarry Smith   ex->relax_type     = 1;
12239e5bc791SBarry Smith   ex->rap_type       = 0;
12249e5bc791SBarry Smith   ex->num_pre_relax  = 1;
12259e5bc791SBarry Smith   ex->num_post_relax = 1;
12263b46a515SGlenn Hammond   ex->max_levels     = 0;
12279e5bc791SBarry Smith 
1228ebc551c0SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_PFMG;
1229ebc551c0SBarry Smith   pc->ops->view            = PCView_PFMG;
1230ebc551c0SBarry Smith   pc->ops->destroy         = PCDestroy_PFMG;
1231f91d8e95SBarry Smith   pc->ops->apply           = PCApply_PFMG;
12329e5bc791SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_PFMG;
123368326731SBarry Smith   pc->ops->setup           = PCSetUp_PFMG;
12342fa5cd67SKarl Rupp 
1235ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ex->hcomm));CHKERRQ(ierr);
1236fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
1237ebc551c0SBarry Smith   PetscFunctionReturn(0);
1238ebc551c0SBarry Smith }
1239d851a50bSGlenn Hammond 
1240325fc9f4SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------------------------*/
1241325fc9f4SBarry Smith 
1242d851a50bSGlenn Hammond /* we know we are working with a HYPRE_SStructMatrix */
1243d851a50bSGlenn Hammond typedef struct {
1244d851a50bSGlenn Hammond   MPI_Comm            hcomm;       /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */
1245d851a50bSGlenn Hammond   HYPRE_SStructSolver ss_solver;
1246d851a50bSGlenn Hammond 
1247d851a50bSGlenn Hammond   /* keep copy of SYSPFMG options used so may view them */
12484ddd07fcSJed Brown   PetscInt its;
1249d851a50bSGlenn Hammond   double   tol;
12504ddd07fcSJed Brown   PetscInt relax_type;
12514ddd07fcSJed Brown   PetscInt num_pre_relax,num_post_relax;
1252d851a50bSGlenn Hammond } PC_SysPFMG;
1253d851a50bSGlenn Hammond 
1254d851a50bSGlenn Hammond #undef __FUNCT__
1255d851a50bSGlenn Hammond #define __FUNCT__ "PCDestroy_SysPFMG"
1256d851a50bSGlenn Hammond PetscErrorCode PCDestroy_SysPFMG(PC pc)
1257d851a50bSGlenn Hammond {
1258d851a50bSGlenn Hammond   PetscErrorCode ierr;
1259d851a50bSGlenn Hammond   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;
1260d851a50bSGlenn Hammond 
1261d851a50bSGlenn Hammond   PetscFunctionBegin;
12622fa5cd67SKarl Rupp   if (ex->ss_solver) PetscStackCallStandard(HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));
1263d851a50bSGlenn Hammond   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1264c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1265d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1266d851a50bSGlenn Hammond }
1267d851a50bSGlenn Hammond 
1268d851a50bSGlenn Hammond static const char *SysPFMGRelaxType[] = {"Weighted-Jacobi","Red/Black-Gauss-Seidel"};
1269d851a50bSGlenn Hammond 
1270d851a50bSGlenn Hammond #undef __FUNCT__
1271d851a50bSGlenn Hammond #define __FUNCT__ "PCView_SysPFMG"
1272d851a50bSGlenn Hammond PetscErrorCode PCView_SysPFMG(PC pc,PetscViewer viewer)
1273d851a50bSGlenn Hammond {
1274d851a50bSGlenn Hammond   PetscErrorCode ierr;
1275ace3abfcSBarry Smith   PetscBool      iascii;
1276d851a50bSGlenn Hammond   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;
1277d851a50bSGlenn Hammond 
1278d851a50bSGlenn Hammond   PetscFunctionBegin;
1279251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1280d851a50bSGlenn Hammond   if (iascii) {
1281d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG preconditioning\n");CHKERRQ(ierr);
1282d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
1283d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
1284d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
1285d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
1286d851a50bSGlenn Hammond   }
1287d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1288d851a50bSGlenn Hammond }
1289d851a50bSGlenn Hammond 
1290d851a50bSGlenn Hammond 
1291d851a50bSGlenn Hammond #undef __FUNCT__
1292d851a50bSGlenn Hammond #define __FUNCT__ "PCSetFromOptions_SysPFMG"
1293d851a50bSGlenn Hammond PetscErrorCode PCSetFromOptions_SysPFMG(PC pc)
1294d851a50bSGlenn Hammond {
1295d851a50bSGlenn Hammond   PetscErrorCode ierr;
1296d851a50bSGlenn Hammond   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;
1297ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
1298d851a50bSGlenn Hammond 
1299d851a50bSGlenn Hammond   PetscFunctionBegin;
1300d851a50bSGlenn Hammond   ierr = PetscOptionsHead("SysPFMG options");CHKERRQ(ierr);
13010298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_syspfmg_print_statistics","Print statistics","HYPRE_SStructSysPFMGSetPrintLevel",flg,&flg,NULL);CHKERRQ(ierr);
1302d851a50bSGlenn Hammond   if (flg) {
1303d851a50bSGlenn Hammond     int level=3;
1304fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructSysPFMGSetPrintLevel,(ex->ss_solver,level));
1305d851a50bSGlenn Hammond   }
13060298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_syspfmg_its","Number of iterations of SysPFMG to use as preconditioner","HYPRE_SStructSysPFMGSetMaxIter",ex->its,&ex->its,NULL);CHKERRQ(ierr);
1307fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(ex->ss_solver,ex->its));
13080298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_syspfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_SStructSysPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,NULL);CHKERRQ(ierr);
1309fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetNumPreRelax,(ex->ss_solver,ex->num_pre_relax));
13100298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_syspfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_SStructSysPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,NULL);CHKERRQ(ierr);
1311fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetNumPostRelax,(ex->ss_solver,ex->num_post_relax));
1312d851a50bSGlenn Hammond 
13130298fd71SBarry Smith   ierr = PetscOptionsReal("-pc_syspfmg_tol","Tolerance of SysPFMG","HYPRE_SStructSysPFMGSetTol",ex->tol,&ex->tol,NULL);CHKERRQ(ierr);
1314fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(ex->ss_solver,ex->tol));
13150298fd71SBarry Smith   ierr = PetscOptionsEList("-pc_syspfmg_relax_type","Relax type for the up and down cycles","HYPRE_SStructSysPFMGSetRelaxType",SysPFMGRelaxType,4,SysPFMGRelaxType[ex->relax_type],&ex->relax_type,NULL);CHKERRQ(ierr);
1316fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetRelaxType,(ex->ss_solver, ex->relax_type));
1317d851a50bSGlenn Hammond   ierr = PetscOptionsTail();CHKERRQ(ierr);
1318d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1319d851a50bSGlenn Hammond }
1320d851a50bSGlenn Hammond 
1321d851a50bSGlenn Hammond #undef __FUNCT__
1322d851a50bSGlenn Hammond #define __FUNCT__ "PCApply_SysPFMG"
1323d851a50bSGlenn Hammond PetscErrorCode PCApply_SysPFMG(PC pc,Vec x,Vec y)
1324d851a50bSGlenn Hammond {
1325d851a50bSGlenn Hammond   PetscErrorCode   ierr;
1326d851a50bSGlenn Hammond   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
1327d851a50bSGlenn Hammond   PetscScalar      *xx,*yy;
13284ddd07fcSJed Brown   PetscInt         ilower[3],iupper[3];
1329d851a50bSGlenn Hammond   Mat_HYPRESStruct *mx     = (Mat_HYPRESStruct*)(pc->pmat->data);
13304ddd07fcSJed Brown   PetscInt         ordering= mx->dofs_order;
13314ddd07fcSJed Brown   PetscInt         nvars   = mx->nvars;
13324ddd07fcSJed Brown   PetscInt         part    = 0;
13334ddd07fcSJed Brown   PetscInt         size;
13344ddd07fcSJed Brown   PetscInt         i;
1335d851a50bSGlenn Hammond 
1336d851a50bSGlenn Hammond   PetscFunctionBegin;
1337aa219208SBarry Smith   ierr       = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1338d851a50bSGlenn Hammond   iupper[0] += ilower[0] - 1;
1339d851a50bSGlenn Hammond   iupper[1] += ilower[1] - 1;
1340d851a50bSGlenn Hammond   iupper[2] += ilower[2] - 1;
1341d851a50bSGlenn Hammond 
1342d851a50bSGlenn Hammond   size = 1;
13432fa5cd67SKarl Rupp   for (i= 0; i< 3; i++) size *= (iupper[i]-ilower[i]+1);
13442fa5cd67SKarl Rupp 
1345d851a50bSGlenn Hammond   /* copy x values over to hypre for variable ordering */
1346d851a50bSGlenn Hammond   if (ordering) {
1347fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
1348d851a50bSGlenn Hammond     ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
13492fa5cd67SKarl Rupp     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,ilower,iupper,i,xx+(size*i)));
1350d851a50bSGlenn Hammond     ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1351fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
1352fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
1353fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1354d851a50bSGlenn Hammond 
1355d851a50bSGlenn Hammond     /* copy solution values back to PETSc */
1356d851a50bSGlenn Hammond     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
13572fa5cd67SKarl Rupp     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,ilower,iupper,i,yy+(size*i)));
1358d851a50bSGlenn Hammond     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1359a65764d7SBarry Smith   } else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
1360d851a50bSGlenn Hammond     PetscScalar *z;
13614ddd07fcSJed Brown     PetscInt    j, k;
1362d851a50bSGlenn Hammond 
1363d851a50bSGlenn Hammond     ierr = PetscMalloc(nvars*size*sizeof(PetscScalar),&z);CHKERRQ(ierr);
1364fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
1365d851a50bSGlenn Hammond     ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1366d851a50bSGlenn Hammond 
1367d851a50bSGlenn Hammond     /* transform nodal to hypre's variable ordering for sys_pfmg */
1368d851a50bSGlenn Hammond     for (i= 0; i< size; i++) {
1369d851a50bSGlenn Hammond       k= i*nvars;
13702fa5cd67SKarl Rupp       for (j= 0; j< nvars; j++) z[j*size+i]= xx[k+j];
1371d851a50bSGlenn Hammond     }
13722fa5cd67SKarl Rupp     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,ilower,iupper,i,z+(size*i)));
1373d851a50bSGlenn Hammond     ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1374fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
1375fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1376d851a50bSGlenn Hammond 
1377d851a50bSGlenn Hammond     /* copy solution values back to PETSc */
1378d851a50bSGlenn Hammond     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
13792fa5cd67SKarl Rupp     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,ilower,iupper,i,z+(size*i)));
1380d851a50bSGlenn Hammond     /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
1381d851a50bSGlenn Hammond     for (i= 0; i< size; i++) {
1382d851a50bSGlenn Hammond       k= i*nvars;
13832fa5cd67SKarl Rupp       for (j= 0; j< nvars; j++) yy[k+j]= z[j*size+i];
1384d851a50bSGlenn Hammond     }
1385d851a50bSGlenn Hammond     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1386d851a50bSGlenn Hammond     ierr = PetscFree(z);CHKERRQ(ierr);
1387d851a50bSGlenn Hammond   }
1388d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1389d851a50bSGlenn Hammond }
1390d851a50bSGlenn Hammond 
1391d851a50bSGlenn Hammond #undef __FUNCT__
1392d851a50bSGlenn Hammond #define __FUNCT__ "PCApplyRichardson_SysPFMG"
1393ace3abfcSBarry Smith 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)
1394d851a50bSGlenn Hammond {
1395d851a50bSGlenn Hammond   PC_SysPFMG     *jac = (PC_SysPFMG*)pc->data;
1396d851a50bSGlenn Hammond   PetscErrorCode ierr;
13974ddd07fcSJed Brown   PetscInt       oits;
1398d851a50bSGlenn Hammond 
1399d851a50bSGlenn Hammond   PetscFunctionBegin;
1400fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,its*jac->its));
1401fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,rtol));
1402d851a50bSGlenn Hammond   ierr = PCApply_SysPFMG(pc,b,y);CHKERRQ(ierr);
1403fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGGetNumIterations,(jac->ss_solver,&oits));
1404d851a50bSGlenn Hammond   *outits = oits;
1405d851a50bSGlenn Hammond   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
1406d851a50bSGlenn Hammond   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
1407fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,jac->tol));
1408fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,jac->its));
1409d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1410d851a50bSGlenn Hammond }
1411d851a50bSGlenn Hammond 
1412d851a50bSGlenn Hammond 
1413d851a50bSGlenn Hammond #undef __FUNCT__
1414d851a50bSGlenn Hammond #define __FUNCT__ "PCSetUp_SysPFMG"
1415d851a50bSGlenn Hammond PetscErrorCode PCSetUp_SysPFMG(PC pc)
1416d851a50bSGlenn Hammond {
1417d851a50bSGlenn Hammond   PetscErrorCode   ierr;
1418d851a50bSGlenn Hammond   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
1419d851a50bSGlenn Hammond   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct*)(pc->pmat->data);
1420ace3abfcSBarry Smith   PetscBool        flg;
1421d851a50bSGlenn Hammond 
1422d851a50bSGlenn Hammond   PetscFunctionBegin;
1423251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESSTRUCT,&flg);CHKERRQ(ierr);
1424ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESSTRUCT with this preconditioner");
1425d851a50bSGlenn Hammond 
1426d851a50bSGlenn Hammond   /* create the hypre sstruct solver object and set its information */
14272fa5cd67SKarl Rupp   if (ex->ss_solver) PetscStackCallStandard(HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));
1428fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
1429d851a50bSGlenn Hammond   ierr = PCSetFromOptions_SysPFMG(pc);CHKERRQ(ierr);
1430fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetZeroGuess,(ex->ss_solver));
1431fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetup,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1432d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1433d851a50bSGlenn Hammond }
1434d851a50bSGlenn Hammond 
1435d851a50bSGlenn Hammond 
1436d851a50bSGlenn Hammond /*MC
1437d851a50bSGlenn Hammond      PCSysPFMG - the hypre SysPFMG multigrid solver
1438d851a50bSGlenn Hammond 
1439d851a50bSGlenn Hammond    Level: advanced
1440d851a50bSGlenn Hammond 
1441d851a50bSGlenn Hammond    Options Database:
1442d851a50bSGlenn Hammond + -pc_syspfmg_its <its> number of iterations of SysPFMG to use as preconditioner
1443d851a50bSGlenn Hammond . -pc_syspfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
1444d851a50bSGlenn Hammond . -pc_syspfmg_num_post_relax <steps> number of smoothing steps after coarse grid
1445d851a50bSGlenn Hammond . -pc_syspfmg_tol <tol> tolerance of SysPFMG
1446d851a50bSGlenn Hammond . -pc_syspfmg_relax_type -relaxation type for the up and down cycles, one of Weighted-Jacobi,Red/Black-Gauss-Seidel
1447d851a50bSGlenn Hammond 
1448d851a50bSGlenn Hammond    Notes:  This is for CELL-centered descretizations
1449d851a50bSGlenn Hammond 
1450f6680f47SSatish Balay            This must be used with the MATHYPRESSTRUCT matrix type.
1451aa219208SBarry Smith            This is less general than in hypre, it supports only one part, and one block per process defined by a PETSc DMDA.
1452d851a50bSGlenn Hammond            Also, only cell-centered variables.
1453d851a50bSGlenn Hammond 
1454d851a50bSGlenn Hammond .seealso:  PCMG, MATHYPRESSTRUCT
1455d851a50bSGlenn Hammond M*/
1456d851a50bSGlenn Hammond 
1457d851a50bSGlenn Hammond #undef __FUNCT__
1458d851a50bSGlenn Hammond #define __FUNCT__ "PCCreate_SysPFMG"
1459*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SysPFMG(PC pc)
1460d851a50bSGlenn Hammond {
1461d851a50bSGlenn Hammond   PetscErrorCode ierr;
1462d851a50bSGlenn Hammond   PC_SysPFMG     *ex;
1463d851a50bSGlenn Hammond 
1464d851a50bSGlenn Hammond   PetscFunctionBegin;
1465d851a50bSGlenn Hammond   ierr     = PetscNew(PC_SysPFMG,&ex);CHKERRQ(ierr); \
1466d851a50bSGlenn Hammond   pc->data = ex;
1467d851a50bSGlenn Hammond 
1468d851a50bSGlenn Hammond   ex->its            = 1;
1469d851a50bSGlenn Hammond   ex->tol            = 1.e-8;
1470d851a50bSGlenn Hammond   ex->relax_type     = 1;
1471d851a50bSGlenn Hammond   ex->num_pre_relax  = 1;
1472d851a50bSGlenn Hammond   ex->num_post_relax = 1;
1473d851a50bSGlenn Hammond 
1474d851a50bSGlenn Hammond   pc->ops->setfromoptions  = PCSetFromOptions_SysPFMG;
1475d851a50bSGlenn Hammond   pc->ops->view            = PCView_SysPFMG;
1476d851a50bSGlenn Hammond   pc->ops->destroy         = PCDestroy_SysPFMG;
1477d851a50bSGlenn Hammond   pc->ops->apply           = PCApply_SysPFMG;
1478d851a50bSGlenn Hammond   pc->ops->applyrichardson = PCApplyRichardson_SysPFMG;
1479d851a50bSGlenn Hammond   pc->ops->setup           = PCSetUp_SysPFMG;
14802fa5cd67SKarl Rupp 
1481ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ex->hcomm));CHKERRQ(ierr);
1482fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
1483d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1484d851a50bSGlenn Hammond }
1485