xref: /petsc/src/ksp/pc/impls/hypre/hypre.c (revision 0ad7597df3318518cf08ebbaa20b6a237bca446d)
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);
12416d9e3a6SLisandro Dalcin     if (bs > 1) {
125fd3f9acdSBarry Smith       PetscStackCallStandard(HYPRE_BoomerAMGSetNumFunctions,(jac->hsolver,bs));
12616d9e3a6SLisandro Dalcin     }
12716d9e3a6SLisandro Dalcin   };
12816d9e3a6SLisandro Dalcin   ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr);
129fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
130fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&bv));
131fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&xv));
132fd3f9acdSBarry Smith   PetscStackCall("HYPRE_SetupXXX",ierr = (*jac->setup)(jac->hsolver,hmat,bv,xv);CHKERRQ(ierr););
13316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
13416d9e3a6SLisandro Dalcin }
13516d9e3a6SLisandro Dalcin 
13616d9e3a6SLisandro Dalcin /*
13716d9e3a6SLisandro Dalcin     Replaces the address where the HYPRE vector points to its data with the address of
13816d9e3a6SLisandro Dalcin   PETSc's data. Saves the old address so it can be reset when we are finished with it.
13916d9e3a6SLisandro Dalcin   Allows use to get the data into a HYPRE vector without the cost of memcopies
14016d9e3a6SLisandro Dalcin */
14116d9e3a6SLisandro Dalcin #define HYPREReplacePointer(b,newvalue,savedvalue) {\
14216d9e3a6SLisandro Dalcin    hypre_ParVector *par_vector   = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));\
14316d9e3a6SLisandro Dalcin    hypre_Vector    *local_vector = hypre_ParVectorLocalVector(par_vector);\
14416d9e3a6SLisandro Dalcin    savedvalue         = local_vector->data;\
145*0ad7597dSKarl Rupp    local_vector->data = newvalue;          \
146*0ad7597dSKarl Rupp    }
14716d9e3a6SLisandro Dalcin 
14816d9e3a6SLisandro Dalcin #undef __FUNCT__
14916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApply_HYPRE"
15016d9e3a6SLisandro Dalcin static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
15116d9e3a6SLisandro Dalcin {
15216d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
15316d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
15416d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
15516d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
15616d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
15716d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
1584ddd07fcSJed Brown   PetscInt           hierr;
15916d9e3a6SLisandro Dalcin 
16016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
16116d9e3a6SLisandro Dalcin   if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);}
16216d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
16316d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
16416d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
16516d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
16616d9e3a6SLisandro Dalcin 
167fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
168fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&jbv));
169fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&jxv));
170fd3f9acdSBarry Smith   PetscStackCall("Hypre solve",hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv);
17165e19b50SBarry Smith                                if (hierr && hierr != HYPRE_ERROR_CONV) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
172fd3f9acdSBarry Smith                                if (hierr) hypre__global_error = 0;);
17316d9e3a6SLisandro Dalcin 
17416d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
17516d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
17616d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
17716d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
17816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
17916d9e3a6SLisandro Dalcin }
18016d9e3a6SLisandro Dalcin 
18116d9e3a6SLisandro Dalcin #undef __FUNCT__
18216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCDestroy_HYPRE"
18316d9e3a6SLisandro Dalcin static PetscErrorCode PCDestroy_HYPRE(PC pc)
18416d9e3a6SLisandro Dalcin {
18516d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
18616d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
18716d9e3a6SLisandro Dalcin 
18816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
189fd3f9acdSBarry Smith   if (jac->ij) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->ij));
190fd3f9acdSBarry Smith   if (jac->b)  PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->b));
191fd3f9acdSBarry Smith   if (jac->x)  PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->x));
192fd3f9acdSBarry Smith   if (jac->destroy) PetscStackCall("HYPRE_DistroyXXX",ierr = (*jac->destroy)(jac->hsolver);CHKERRQ(ierr););
193503cfb0cSBarry Smith   ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr);
19416d9e3a6SLisandro Dalcin   if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);}
195c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
19616d9e3a6SLisandro Dalcin 
19716d9e3a6SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr);
19816d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);CHKERRQ(ierr);
19916d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);CHKERRQ(ierr);
20016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
20116d9e3a6SLisandro Dalcin }
20216d9e3a6SLisandro Dalcin 
20316d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
20416d9e3a6SLisandro Dalcin #undef __FUNCT__
20516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut"
20616d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc)
20716d9e3a6SLisandro Dalcin {
20816d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
20916d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
210ace3abfcSBarry Smith   PetscBool      flag;
21116d9e3a6SLisandro Dalcin 
21216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
21316d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr);
21416d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr);
215fd3f9acdSBarry Smith   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetMaxIter,(jac->hsolver,jac->maxiter));
21616d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr);
217fd3f9acdSBarry Smith   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetDropTolerance,(jac->hsolver,jac->tol));
21816d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr);
219fd3f9acdSBarry Smith   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetFactorRowSize,(jac->hsolver,jac->factorrowsize));
22016d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
22116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
22216d9e3a6SLisandro Dalcin }
22316d9e3a6SLisandro Dalcin 
22416d9e3a6SLisandro Dalcin #undef __FUNCT__
22516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Pilut"
22616d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
22716d9e3a6SLisandro Dalcin {
22816d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
22916d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
230ace3abfcSBarry Smith   PetscBool      iascii;
23116d9e3a6SLisandro Dalcin 
23216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
233251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
23416d9e3a6SLisandro Dalcin   if (iascii) {
23516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");CHKERRQ(ierr);
23616d9e3a6SLisandro Dalcin     if (jac->maxiter != PETSC_DEFAULT) {
23716d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr);
23816d9e3a6SLisandro Dalcin     } else {
23916d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr);
24016d9e3a6SLisandro Dalcin     }
24116d9e3a6SLisandro Dalcin     if (jac->tol != PETSC_DEFAULT) {
24216d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr);
24316d9e3a6SLisandro Dalcin     } else {
24416d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr);
24516d9e3a6SLisandro Dalcin     }
24616d9e3a6SLisandro Dalcin     if (jac->factorrowsize != PETSC_DEFAULT) {
24716d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr);
24816d9e3a6SLisandro Dalcin     } else {
24916d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");CHKERRQ(ierr);
25016d9e3a6SLisandro Dalcin     }
25116d9e3a6SLisandro Dalcin   }
25216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
25316d9e3a6SLisandro Dalcin }
25416d9e3a6SLisandro Dalcin 
25516d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
25616d9e3a6SLisandro Dalcin #undef __FUNCT__
25716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid"
25816d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc)
25916d9e3a6SLisandro Dalcin {
26016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
26116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
262ace3abfcSBarry Smith   PetscBool      flag;
263390e7148SBarry Smith   char           *args[8],levels[16];
264390e7148SBarry Smith   PetscInt       cnt = 0;
26516d9e3a6SLisandro Dalcin 
26616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
26716d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr);
26816d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr);
26916d9e3a6SLisandro Dalcin   if (flag) {
27065e19b50SBarry Smith     if (jac->levels < 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels);
2718caf3d72SBarry Smith     ierr = PetscSNPrintf(levels,sizeof(levels),"%D",jac->levels);CHKERRQ(ierr);
272390e7148SBarry Smith     args[cnt++] = (char*)"-level"; args[cnt++] = levels;
27316d9e3a6SLisandro Dalcin   }
274acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);CHKERRQ(ierr);
27516d9e3a6SLisandro Dalcin   if (jac->bjilu) {
276390e7148SBarry Smith     args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1";
27716d9e3a6SLisandro Dalcin   }
27816d9e3a6SLisandro Dalcin 
279acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);CHKERRQ(ierr);
28016d9e3a6SLisandro Dalcin   if (jac->printstatistics) {
281390e7148SBarry Smith     args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1";
282390e7148SBarry Smith     args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1";
28316d9e3a6SLisandro Dalcin   }
28416d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
285fd3f9acdSBarry Smith   if (cnt) PetscStackCallStandard(HYPRE_EuclidSetParams,(jac->hsolver,cnt,args));
28616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
28716d9e3a6SLisandro Dalcin }
28816d9e3a6SLisandro Dalcin 
28916d9e3a6SLisandro Dalcin #undef __FUNCT__
29016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Euclid"
29116d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
29216d9e3a6SLisandro Dalcin {
29316d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
29416d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
295ace3abfcSBarry Smith   PetscBool      iascii;
29616d9e3a6SLisandro Dalcin 
29716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
298251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
29916d9e3a6SLisandro Dalcin   if (iascii) {
30016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n");CHKERRQ(ierr);
30116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr);
30216d9e3a6SLisandro Dalcin     if (jac->bjilu) {
30316d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr);
30416d9e3a6SLisandro Dalcin     }
30516d9e3a6SLisandro Dalcin   }
30616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
30716d9e3a6SLisandro Dalcin }
30816d9e3a6SLisandro Dalcin 
30916d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
31016d9e3a6SLisandro Dalcin 
31116d9e3a6SLisandro Dalcin #undef __FUNCT__
31216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG"
31316d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
31416d9e3a6SLisandro Dalcin {
31516d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
31616d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
31716d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
31816d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
31916d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
32016d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
3214ddd07fcSJed Brown   PetscInt           hierr;
32216d9e3a6SLisandro Dalcin 
32316d9e3a6SLisandro Dalcin   PetscFunctionBegin;
32416d9e3a6SLisandro Dalcin   ierr = VecSet(x,0.0);CHKERRQ(ierr);
32516d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
32616d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
32716d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
32816d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
32916d9e3a6SLisandro Dalcin 
330fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
331fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&jbv));
332fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&jxv));
33316d9e3a6SLisandro Dalcin 
33416d9e3a6SLisandro Dalcin   hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv);
33516d9e3a6SLisandro Dalcin   /* error code of 1 in BoomerAMG merely means convergence not achieved */
336e32f2f54SBarry Smith   if (hierr && (hierr != 1)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
33716d9e3a6SLisandro Dalcin   if (hierr) hypre__global_error = 0;
33816d9e3a6SLisandro Dalcin 
33916d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
34016d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
34116d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
34216d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
34316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
34416d9e3a6SLisandro Dalcin }
34516d9e3a6SLisandro Dalcin 
346a669f990SJed Brown /* static array length */
347a669f990SJed Brown #define ALEN(a) (sizeof(a)/sizeof((a)[0]))
348a669f990SJed Brown 
34916d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
3500f1074feSSatish Balay static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"};
35116d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
35265de4495SJed Brown /* The following corresponds to HYPRE_BoomerAMGSetRelaxType which has many missing numbers in the enum */
35365de4495SJed Brown static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","seqboundary-Gauss-Seidel","SOR/Jacobi","backward-SOR/Jacobi",
35465de4495SJed Brown                                                   "" /* [5] hybrid chaotic Gauss-Seidel (works only with OpenMP) */,"symmetric-SOR/Jacobi",
35565de4495SJed Brown                                                   "" /* 7 */,"l1scaled-SOR/Jacobi","Gaussian-elimination",
35665de4495SJed Brown                                                   "" /* 10 */, "" /* 11 */, "" /* 12 */, "" /* 13 */, "" /* 14 */,
35765de4495SJed Brown                                                   "CG" /* non-stationary */,"Chebyshev","FCF-Jacobi","l1scaled-Jacobi"};
3580f1074feSSatish Balay static const char *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i",
3590f1074feSSatish Balay                                                   "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"};
36016d9e3a6SLisandro Dalcin #undef __FUNCT__
36116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG"
36216d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
36316d9e3a6SLisandro Dalcin {
36416d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
36516d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
3664ddd07fcSJed Brown   PetscInt       n,indx,level;
367ace3abfcSBarry Smith   PetscBool      flg, tmp_truth;
36816d9e3a6SLisandro Dalcin   double         tmpdbl, twodbl[2];
36916d9e3a6SLisandro Dalcin 
37016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
37116d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr);
3724336a9eeSBarry Smith   ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType+1,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr);
37316d9e3a6SLisandro Dalcin   if (flg) {
3744336a9eeSBarry Smith     jac->cycletype = indx+1;
375fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
37616d9e3a6SLisandro Dalcin   }
37716d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr);
37816d9e3a6SLisandro Dalcin   if (flg) {
37965e19b50SBarry Smith     if (jac->maxlevels < 2) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
380fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
38116d9e3a6SLisandro Dalcin   }
38216d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr);
38316d9e3a6SLisandro Dalcin   if (flg) {
38465e19b50SBarry Smith     if (jac->maxiter < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
385fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
38616d9e3a6SLisandro Dalcin   }
3870f1074feSSatish 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);
38816d9e3a6SLisandro Dalcin   if (flg) {
38965e19b50SBarry Smith     if (jac->tol < 0.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be greater than or equal to zero",jac->tol);
390fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
39116d9e3a6SLisandro Dalcin   }
39216d9e3a6SLisandro Dalcin 
3930f1074feSSatish Balay   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr);
39416d9e3a6SLisandro Dalcin   if (flg) {
39565e19b50SBarry Smith     if (jac->truncfactor < 0.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor);
396fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
39716d9e3a6SLisandro Dalcin   }
39816d9e3a6SLisandro Dalcin 
3990f1074feSSatish 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);
4000f1074feSSatish Balay   if (flg) {
40165e19b50SBarry Smith     if (jac->pmax < 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"P_max %G must be greater than or equal to zero",jac->pmax);
402fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
4030f1074feSSatish Balay   }
4040f1074feSSatish Balay 
4050f1074feSSatish Balay  ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);CHKERRQ(ierr);
4060f1074feSSatish Balay   if (flg) {
40765e19b50SBarry Smith      if (jac->agg_nl < 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %G must be greater than or equal to zero",jac->agg_nl);
4080f1074feSSatish Balay 
409fd3f9acdSBarry Smith      PetscStackCallStandard(HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
4100f1074feSSatish Balay   }
4110f1074feSSatish Balay 
4120f1074feSSatish Balay 
4130f1074feSSatish 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);
4140f1074feSSatish Balay   if (flg) {
41565e19b50SBarry Smith      if (jac->agg_num_paths < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of paths %G must be greater than or equal to 1",jac->agg_num_paths);
4160f1074feSSatish Balay 
417fd3f9acdSBarry Smith      PetscStackCallStandard(HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
4180f1074feSSatish Balay   }
4190f1074feSSatish Balay 
4200f1074feSSatish Balay 
42116d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr);
42216d9e3a6SLisandro Dalcin   if (flg) {
42365e19b50SBarry Smith     if (jac->strongthreshold < 0.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold);
424fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
42516d9e3a6SLisandro Dalcin   }
42616d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr);
42716d9e3a6SLisandro Dalcin   if (flg) {
42865e19b50SBarry Smith     if (jac->maxrowsum < 0.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum);
42965e19b50SBarry Smith     if (jac->maxrowsum > 1.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum);
430fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
43116d9e3a6SLisandro Dalcin   }
43216d9e3a6SLisandro Dalcin 
43316d9e3a6SLisandro Dalcin   /* Grid sweeps */
4340f1074feSSatish 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);
43516d9e3a6SLisandro Dalcin   if (flg) {
436fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver,indx));
43716d9e3a6SLisandro Dalcin     /* modify the jac structure so we can view the updated options with PC_View */
43816d9e3a6SLisandro Dalcin     jac->gridsweeps[0] = indx;
4390f1074feSSatish Balay     jac->gridsweeps[1] = indx;
4400f1074feSSatish Balay     /*defaults coarse to 1 */
4410f1074feSSatish Balay     jac->gridsweeps[2] = 1;
44216d9e3a6SLisandro Dalcin   }
4430f1074feSSatish Balay 
4440f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx ,&flg);CHKERRQ(ierr);
44516d9e3a6SLisandro Dalcin   if (flg) {
446fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 1));
4470f1074feSSatish Balay     jac->gridsweeps[0] = indx;
44816d9e3a6SLisandro Dalcin   }
44916d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr);
45016d9e3a6SLisandro Dalcin   if (flg) {
451fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 2));
4520f1074feSSatish Balay     jac->gridsweeps[1] = indx;
45316d9e3a6SLisandro Dalcin   }
4540f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr);
45516d9e3a6SLisandro Dalcin   if (flg) {
456fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 3));
4570f1074feSSatish Balay     jac->gridsweeps[2] = indx;
45816d9e3a6SLisandro Dalcin   }
45916d9e3a6SLisandro Dalcin 
46016d9e3a6SLisandro Dalcin   /* Relax type */
461a669f990SJed 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);
46216d9e3a6SLisandro Dalcin   if (flg) {
4630f1074feSSatish Balay     jac->relaxtype[0] = jac->relaxtype[1]  = indx;
464fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, indx));
4650f1074feSSatish Balay     /* by default, coarse type set to 9 */
4660f1074feSSatish Balay     jac->relaxtype[2] = 9;
4670f1074feSSatish Balay 
46816d9e3a6SLisandro Dalcin   }
469a669f990SJed 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);
47016d9e3a6SLisandro Dalcin   if (flg) {
47116d9e3a6SLisandro Dalcin     jac->relaxtype[0] = indx;
472fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 1));
47316d9e3a6SLisandro Dalcin   }
474a669f990SJed 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);
47516d9e3a6SLisandro Dalcin   if (flg) {
4760f1074feSSatish Balay     jac->relaxtype[1] = indx;
477fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 2));
47816d9e3a6SLisandro Dalcin   }
479a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr);
48016d9e3a6SLisandro Dalcin   if (flg) {
4810f1074feSSatish Balay     jac->relaxtype[2] = indx;
482fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 3));
48316d9e3a6SLisandro Dalcin   }
48416d9e3a6SLisandro Dalcin 
48516d9e3a6SLisandro Dalcin   /* Relaxation Weight */
48616d9e3a6SLisandro 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);
48716d9e3a6SLisandro Dalcin   if (flg) {
488fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxWt,(jac->hsolver,tmpdbl));
48916d9e3a6SLisandro Dalcin     jac->relaxweight = tmpdbl;
49016d9e3a6SLisandro Dalcin   }
49116d9e3a6SLisandro Dalcin 
49216d9e3a6SLisandro Dalcin   n=2;
49316d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
49416d9e3a6SLisandro 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);
49516d9e3a6SLisandro Dalcin   if (flg) {
49616d9e3a6SLisandro Dalcin     if (n == 2) {
49716d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
498fd3f9acdSBarry Smith       PetscStackCallStandard(HYPRE_BoomerAMGSetLevelRelaxWt,(jac->hsolver,twodbl[0],indx));
49916d9e3a6SLisandro Dalcin     } else {
50065e19b50SBarry Smith       SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relax weight level: you must provide 2 values separated by a comma (and no space), you provided %d",n);
50116d9e3a6SLisandro Dalcin     }
50216d9e3a6SLisandro Dalcin   }
50316d9e3a6SLisandro Dalcin 
50416d9e3a6SLisandro Dalcin   /* Outer relaxation Weight */
50516d9e3a6SLisandro 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);
50616d9e3a6SLisandro Dalcin   if (flg) {
507fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetOuterWt,(jac->hsolver, tmpdbl));
50816d9e3a6SLisandro Dalcin     jac->outerrelaxweight = tmpdbl;
50916d9e3a6SLisandro Dalcin   }
51016d9e3a6SLisandro Dalcin 
51116d9e3a6SLisandro Dalcin   n=2;
51216d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
51316d9e3a6SLisandro 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);
51416d9e3a6SLisandro Dalcin   if (flg) {
51516d9e3a6SLisandro Dalcin     if (n == 2) {
51616d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
517fd3f9acdSBarry Smith       PetscStackCallStandard(HYPRE_BoomerAMGSetLevelOuterWt,(jac->hsolver, twodbl[0], indx));
51816d9e3a6SLisandro Dalcin     } else {
51965e19b50SBarry Smith       SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relax weight outer level: You must provide 2 values separated by a comma (and no space), you provided %d",n);
52016d9e3a6SLisandro Dalcin     }
52116d9e3a6SLisandro Dalcin   }
52216d9e3a6SLisandro Dalcin 
52316d9e3a6SLisandro Dalcin   /* the Relax Order */
524acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
52516d9e3a6SLisandro Dalcin 
52616d9e3a6SLisandro Dalcin   if (flg) {
52716d9e3a6SLisandro Dalcin     jac->relaxorder = 0;
528fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
52916d9e3a6SLisandro Dalcin   }
530a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,ALEN(HYPREBoomerAMGMeasureType),HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr);
53116d9e3a6SLisandro Dalcin   if (flg) {
53216d9e3a6SLisandro Dalcin     jac->measuretype = indx;
533fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
53416d9e3a6SLisandro Dalcin   }
5350f1074feSSatish Balay   /* update list length 3/07 */
536a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,ALEN(HYPREBoomerAMGCoarsenType),HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr);
53716d9e3a6SLisandro Dalcin   if (flg) {
53816d9e3a6SLisandro Dalcin     jac->coarsentype = indx;
539fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
54016d9e3a6SLisandro Dalcin   }
5410f1074feSSatish Balay 
5420f1074feSSatish Balay   /* new 3/07 */
543a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,ALEN(HYPREBoomerAMGInterpType),HYPREBoomerAMGInterpType[0],&indx,&flg);CHKERRQ(ierr);
5440f1074feSSatish Balay   if (flg) {
5450f1074feSSatish Balay     jac->interptype = indx;
546fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
5470f1074feSSatish Balay   }
5480f1074feSSatish Balay 
549b96a4a96SBarry Smith   ierr = PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);CHKERRQ(ierr);
55016d9e3a6SLisandro Dalcin   if (flg) {
551b96a4a96SBarry Smith     level = 3;
55216d9e3a6SLisandro Dalcin     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
553b96a4a96SBarry Smith     jac->printstatistics = PETSC_TRUE;
554fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetPrintLevel,(jac->hsolver,level));
5552ae77aedSBarry Smith   }
5562ae77aedSBarry Smith 
557b96a4a96SBarry Smith   ierr = PetscOptionsName("-pc_hypre_boomeramg_print_debug","Print debug information","None",&flg);CHKERRQ(ierr);
5582ae77aedSBarry Smith   if (flg) {
559b96a4a96SBarry Smith     level = 3;
5602ae77aedSBarry Smith     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
561b96a4a96SBarry Smith     jac->printstatistics = PETSC_TRUE;
562fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetDebugFlag,(jac->hsolver,level));
56316d9e3a6SLisandro Dalcin   }
5648f87f92bSBarry Smith 
565acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5668f87f92bSBarry Smith   if (flg && tmp_truth) {
5678f87f92bSBarry Smith     jac->nodal_coarsen = 1;
568fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetNodal,(jac->hsolver,1));
5698f87f92bSBarry Smith   }
5708f87f92bSBarry Smith 
571acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5728f87f92bSBarry Smith   if (flg && tmp_truth) {
5738f87f92bSBarry Smith     PetscInt tmp_int;
5748f87f92bSBarry Smith     ierr = PetscOptionsInt("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);CHKERRQ(ierr);
5758f87f92bSBarry Smith     if (flg) jac->nodal_relax_levels = tmp_int;
576fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothType,(jac->hsolver,6));
577fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetDomainType,(jac->hsolver,1));
578fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetOverlap,(jac->hsolver,0));
579fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothNumLevels,(jac->hsolver,jac->nodal_relax_levels));
5808f87f92bSBarry Smith   }
5818f87f92bSBarry Smith 
58216d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
58316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
58416d9e3a6SLisandro Dalcin }
58516d9e3a6SLisandro Dalcin 
58616d9e3a6SLisandro Dalcin #undef __FUNCT__
58716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG"
588ace3abfcSBarry 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)
58916d9e3a6SLisandro Dalcin {
59016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
59116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
5924ddd07fcSJed Brown   PetscInt       oits;
59316d9e3a6SLisandro Dalcin 
59416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
595fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,its*jac->maxiter));
596fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,rtol));
59716d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_TRUE;
59816d9e3a6SLisandro Dalcin   ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr);
59916d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_FALSE;
600fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGGetNumIterations,(jac->hsolver,&oits));
6014d0a8057SBarry Smith   *outits = oits;
6024d0a8057SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
6034d0a8057SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
604fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
605fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
60616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
60716d9e3a6SLisandro Dalcin }
60816d9e3a6SLisandro Dalcin 
60916d9e3a6SLisandro Dalcin 
61016d9e3a6SLisandro Dalcin #undef __FUNCT__
61116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_BoomerAMG"
61216d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
61316d9e3a6SLisandro Dalcin {
61416d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
61516d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
616ace3abfcSBarry Smith   PetscBool      iascii;
61716d9e3a6SLisandro Dalcin 
61816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
619251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
62016d9e3a6SLisandro Dalcin   if (iascii) {
62116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr);
62216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr);
62316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr);
62416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr);
62516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr);
62616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr);
6270f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation truncation factor %G\n",jac->truncfactor);CHKERRQ(ierr);
6280f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);CHKERRQ(ierr);
6290f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);CHKERRQ(ierr);
6300f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);CHKERRQ(ierr);
6310f1074feSSatish Balay 
63216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr);
63316d9e3a6SLisandro Dalcin 
6340f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps down         %d\n",jac->gridsweeps[0]);CHKERRQ(ierr);
6350f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps up           %d\n",jac->gridsweeps[1]);CHKERRQ(ierr);
6360f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on coarse    %d\n",jac->gridsweeps[2]);CHKERRQ(ierr);
63716d9e3a6SLisandro Dalcin 
6380f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr);
6390f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr);
6400f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr);
64116d9e3a6SLisandro Dalcin 
64216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax weight  (all)      %G\n",jac->relaxweight);CHKERRQ(ierr);
64316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr);
64416d9e3a6SLisandro Dalcin 
64516d9e3a6SLisandro Dalcin     if (jac->relaxorder) {
64616d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr);
64716d9e3a6SLisandro Dalcin     } else {
64816d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr);
64916d9e3a6SLisandro Dalcin     }
65016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr);
65116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr);
6520f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation type  %s\n",HYPREBoomerAMGInterpType[jac->interptype]);CHKERRQ(ierr);
6538f87f92bSBarry Smith     if (jac->nodal_coarsen) {
6548f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal())\n");CHKERRQ(ierr);
6558f87f92bSBarry Smith     }
6568f87f92bSBarry Smith     if (jac->nodal_relax) {
6578f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);CHKERRQ(ierr);
6588f87f92bSBarry Smith     }
65916d9e3a6SLisandro Dalcin   }
66016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
66116d9e3a6SLisandro Dalcin }
66216d9e3a6SLisandro Dalcin 
66316d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
66416d9e3a6SLisandro Dalcin #undef __FUNCT__
66516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails"
66616d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc)
66716d9e3a6SLisandro Dalcin {
66816d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
66916d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
6704ddd07fcSJed Brown   PetscInt       indx;
671ace3abfcSBarry Smith   PetscBool      flag;
67216d9e3a6SLisandro Dalcin   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};
67316d9e3a6SLisandro Dalcin 
67416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
67516d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr);
67616d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr);
67716d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr);
67816d9e3a6SLisandro Dalcin   if (flag) {
679fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels));
68016d9e3a6SLisandro Dalcin   }
68116d9e3a6SLisandro Dalcin 
68216d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr);
68316d9e3a6SLisandro Dalcin   if (flag) {
684fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter));
68516d9e3a6SLisandro Dalcin   }
68616d9e3a6SLisandro Dalcin 
68716d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr);
68816d9e3a6SLisandro Dalcin   if (flag) {
689fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal));
69016d9e3a6SLisandro Dalcin   }
69116d9e3a6SLisandro Dalcin 
692acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_parasails_logging","Print info to screen","None",(PetscBool)jac->logging,(PetscBool*)&jac->logging,&flag);CHKERRQ(ierr);
69316d9e3a6SLisandro Dalcin   if (flag) {
694fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging));
69516d9e3a6SLisandro Dalcin   }
69616d9e3a6SLisandro Dalcin 
697acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscBool)jac->ruse,(PetscBool*)&jac->ruse,&flag);CHKERRQ(ierr);
69816d9e3a6SLisandro Dalcin   if (flag) {
699fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse));
70016d9e3a6SLisandro Dalcin   }
70116d9e3a6SLisandro Dalcin 
702a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,ALEN(symtlist),symtlist[0],&indx,&flag);CHKERRQ(ierr);
70316d9e3a6SLisandro Dalcin   if (flag) {
70416d9e3a6SLisandro Dalcin     jac->symt = indx;
705fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
70616d9e3a6SLisandro Dalcin   }
70716d9e3a6SLisandro Dalcin 
70816d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
70916d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
71016d9e3a6SLisandro Dalcin }
71116d9e3a6SLisandro Dalcin 
71216d9e3a6SLisandro Dalcin #undef __FUNCT__
71316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_ParaSails"
71416d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
71516d9e3a6SLisandro Dalcin {
71616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
71716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
718ace3abfcSBarry Smith   PetscBool      iascii;
71916d9e3a6SLisandro Dalcin   const char     *symt = 0;;
72016d9e3a6SLisandro Dalcin 
72116d9e3a6SLisandro Dalcin   PetscFunctionBegin;
722251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
72316d9e3a6SLisandro Dalcin   if (iascii) {
72416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");CHKERRQ(ierr);
72516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr);
72616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr);
72716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr);
72816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr);
729ace3abfcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscBools[jac->ruse]);CHKERRQ(ierr);
730ace3abfcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscBools[jac->logging]);CHKERRQ(ierr);
73116d9e3a6SLisandro Dalcin     if (!jac->symt) {
73216d9e3a6SLisandro Dalcin       symt = "nonsymmetric matrix and preconditioner";
73316d9e3a6SLisandro Dalcin     } else if (jac->symt == 1) {
73416d9e3a6SLisandro Dalcin       symt = "SPD matrix and preconditioner";
73516d9e3a6SLisandro Dalcin     } else if (jac->symt == 2) {
73616d9e3a6SLisandro Dalcin       symt = "nonsymmetric matrix but SPD preconditioner";
73716d9e3a6SLisandro Dalcin     } else {
73865e19b50SBarry Smith       SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
73916d9e3a6SLisandro Dalcin     }
74016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr);
74116d9e3a6SLisandro Dalcin   }
74216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
74316d9e3a6SLisandro Dalcin }
74416d9e3a6SLisandro Dalcin /* ---------------------------------------------------------------------------------*/
74516d9e3a6SLisandro Dalcin 
74616d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
74716d9e3a6SLisandro Dalcin #undef __FUNCT__
74816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType_HYPRE"
7497087cfbeSBarry Smith PetscErrorCode  PCHYPREGetType_HYPRE(PC pc,const char *name[])
75016d9e3a6SLisandro Dalcin {
75116d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
75216d9e3a6SLisandro Dalcin 
75316d9e3a6SLisandro Dalcin   PetscFunctionBegin;
75416d9e3a6SLisandro Dalcin   *name = jac->hypre_type;
75516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
75616d9e3a6SLisandro Dalcin }
75716d9e3a6SLisandro Dalcin EXTERN_C_END
75816d9e3a6SLisandro Dalcin 
75916d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
76016d9e3a6SLisandro Dalcin #undef __FUNCT__
76116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType_HYPRE"
7627087cfbeSBarry Smith PetscErrorCode  PCHYPRESetType_HYPRE(PC pc,const char name[])
76316d9e3a6SLisandro Dalcin {
76416d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
76516d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
766ace3abfcSBarry Smith   PetscBool      flag;
76716d9e3a6SLisandro Dalcin 
76816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
76916d9e3a6SLisandro Dalcin   if (jac->hypre_type) {
77016d9e3a6SLisandro Dalcin     ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr);
771e7e72b3dSBarry Smith     if (!flag) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
77216d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
77316d9e3a6SLisandro Dalcin   } else {
77416d9e3a6SLisandro Dalcin     ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr);
77516d9e3a6SLisandro Dalcin   }
77616d9e3a6SLisandro Dalcin 
77716d9e3a6SLisandro Dalcin   jac->maxiter            = PETSC_DEFAULT;
77816d9e3a6SLisandro Dalcin   jac->tol                = PETSC_DEFAULT;
77916d9e3a6SLisandro Dalcin   jac->printstatistics    = PetscLogPrintInfo;
78016d9e3a6SLisandro Dalcin 
78116d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr);
78216d9e3a6SLisandro Dalcin   if (flag) {
783fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParCSRPilutCreate,(jac->comm_hypre,&jac->hsolver));
78416d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
78516d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Pilut;
78616d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParCSRPilutDestroy;
78716d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParCSRPilutSetup;
78816d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParCSRPilutSolve;
78916d9e3a6SLisandro Dalcin     jac->factorrowsize      = PETSC_DEFAULT;
79016d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
79116d9e3a6SLisandro Dalcin   }
79216d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr);
79316d9e3a6SLisandro Dalcin   if (flag) {
794fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsCreate,(jac->comm_hypre,&jac->hsolver));
79516d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
79616d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_ParaSails;
79716d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParaSailsDestroy;
79816d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParaSailsSetup;
79916d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParaSailsSolve;
80016d9e3a6SLisandro Dalcin     /* initialize */
80116d9e3a6SLisandro Dalcin     jac->nlevels     = 1;
80216d9e3a6SLisandro Dalcin     jac->threshhold  = .1;
80316d9e3a6SLisandro Dalcin     jac->filter      = .1;
80416d9e3a6SLisandro Dalcin     jac->loadbal     = 0;
80516d9e3a6SLisandro Dalcin     if (PetscLogPrintInfo) {
80616d9e3a6SLisandro Dalcin       jac->logging     = (int) PETSC_TRUE;
80716d9e3a6SLisandro Dalcin     } else {
80816d9e3a6SLisandro Dalcin       jac->logging     = (int) PETSC_FALSE;
80916d9e3a6SLisandro Dalcin     }
81016d9e3a6SLisandro Dalcin     jac->ruse = (int) PETSC_FALSE;
81116d9e3a6SLisandro Dalcin     jac->symt   = 0;
812fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels));
813fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter));
814fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal));
815fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging));
816fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse));
817fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
81816d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
81916d9e3a6SLisandro Dalcin   }
82016d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr);
82116d9e3a6SLisandro Dalcin   if (flag) {
82216d9e3a6SLisandro Dalcin     ierr                    = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver);
82316d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
82416d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Euclid;
82516d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_EuclidDestroy;
82616d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_EuclidSetup;
82716d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_EuclidSolve;
82816d9e3a6SLisandro Dalcin     /* initialization */
82916d9e3a6SLisandro Dalcin     jac->bjilu              = PETSC_FALSE;
83016d9e3a6SLisandro Dalcin     jac->levels             = 1;
83116d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
83216d9e3a6SLisandro Dalcin   }
83316d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr);
83416d9e3a6SLisandro Dalcin   if (flag) {
83516d9e3a6SLisandro Dalcin     ierr                     = HYPRE_BoomerAMGCreate(&jac->hsolver);
83616d9e3a6SLisandro Dalcin     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
83716d9e3a6SLisandro Dalcin     pc->ops->view            = PCView_HYPRE_BoomerAMG;
83816d9e3a6SLisandro Dalcin     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
83916d9e3a6SLisandro Dalcin     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
84016d9e3a6SLisandro Dalcin     jac->destroy             = HYPRE_BoomerAMGDestroy;
84116d9e3a6SLisandro Dalcin     jac->setup               = HYPRE_BoomerAMGSetup;
84216d9e3a6SLisandro Dalcin     jac->solve               = HYPRE_BoomerAMGSolve;
84316d9e3a6SLisandro Dalcin     jac->applyrichardson     = PETSC_FALSE;
84416d9e3a6SLisandro Dalcin     /* these defaults match the hypre defaults */
84516d9e3a6SLisandro Dalcin     jac->cycletype        = 1;
84616d9e3a6SLisandro Dalcin     jac->maxlevels        = 25;
84716d9e3a6SLisandro Dalcin     jac->maxiter          = 1;
8488f87f92bSBarry Smith     jac->tol              = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */
84916d9e3a6SLisandro Dalcin     jac->truncfactor      = 0.0;
85016d9e3a6SLisandro Dalcin     jac->strongthreshold  = .25;
85116d9e3a6SLisandro Dalcin     jac->maxrowsum        = .9;
85216d9e3a6SLisandro Dalcin     jac->coarsentype      = 6;
85316d9e3a6SLisandro Dalcin     jac->measuretype      = 0;
8540f1074feSSatish Balay     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = 1;
8558f87f92bSBarry Smith     jac->relaxtype[0]     = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */
8560f1074feSSatish Balay     jac->relaxtype[2]     = 9; /*G.E. */
85716d9e3a6SLisandro Dalcin     jac->relaxweight      = 1.0;
85816d9e3a6SLisandro Dalcin     jac->outerrelaxweight = 1.0;
85916d9e3a6SLisandro Dalcin     jac->relaxorder       = 1;
8600f1074feSSatish Balay     jac->interptype       = 0;
8610f1074feSSatish Balay     jac->agg_nl           = 0;
8620f1074feSSatish Balay     jac->pmax             = 0;
8630f1074feSSatish Balay     jac->truncfactor      = 0.0;
8640f1074feSSatish Balay     jac->agg_num_paths    = 1;
8658f87f92bSBarry Smith 
8668f87f92bSBarry Smith     jac->nodal_coarsen    = 0;
8678f87f92bSBarry Smith     jac->nodal_relax      = PETSC_FALSE;
8688f87f92bSBarry Smith     jac->nodal_relax_levels = 1;
869fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
870fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
871fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
872fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
873fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
874fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
875fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
876fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
877fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
878fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
879fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
880fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
881fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
882fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
883fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, jac->relaxtype[0]));  /*defaults coarse to 9*/
884fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver, jac->gridsweeps[0])); /*defaults coarse to 1 */
88516d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
88616d9e3a6SLisandro Dalcin   }
887503cfb0cSBarry Smith   ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr);
88816d9e3a6SLisandro Dalcin   jac->hypre_type = PETSC_NULL;
88965e19b50SBarry Smith   SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
89016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
89116d9e3a6SLisandro Dalcin }
89216d9e3a6SLisandro Dalcin EXTERN_C_END
89316d9e3a6SLisandro Dalcin 
89416d9e3a6SLisandro Dalcin /*
89516d9e3a6SLisandro Dalcin     It only gets here if the HYPRE type has not been set before the call to
89616d9e3a6SLisandro Dalcin    ...SetFromOptions() which actually is most of the time
89716d9e3a6SLisandro Dalcin */
89816d9e3a6SLisandro Dalcin #undef __FUNCT__
89916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE"
90016d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE(PC pc)
90116d9e3a6SLisandro Dalcin {
90216d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
9034ddd07fcSJed Brown   PetscInt       indx;
90416d9e3a6SLisandro Dalcin   const char     *type[] = {"pilut","parasails","boomeramg","euclid"};
905ace3abfcSBarry Smith   PetscBool      flg;
90616d9e3a6SLisandro Dalcin 
90716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
90816d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr);
90902a17cd4SBarry Smith   ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"boomeramg",&indx,&flg);CHKERRQ(ierr);
91016d9e3a6SLisandro Dalcin   if (flg) {
91116d9e3a6SLisandro Dalcin     ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr);
91202a17cd4SBarry Smith   } else {
91302a17cd4SBarry Smith     ierr = PCHYPRESetType_HYPRE(pc,"boomeramg");CHKERRQ(ierr);
91416d9e3a6SLisandro Dalcin   }
91516d9e3a6SLisandro Dalcin   if (pc->ops->setfromoptions) {
91616d9e3a6SLisandro Dalcin     ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr);
91716d9e3a6SLisandro Dalcin   }
91816d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
91916d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
92016d9e3a6SLisandro Dalcin }
92116d9e3a6SLisandro Dalcin 
92216d9e3a6SLisandro Dalcin #undef __FUNCT__
92316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType"
92416d9e3a6SLisandro Dalcin /*@C
92516d9e3a6SLisandro Dalcin      PCHYPRESetType - Sets which hypre preconditioner you wish to use
92616d9e3a6SLisandro Dalcin 
92716d9e3a6SLisandro Dalcin    Input Parameters:
92816d9e3a6SLisandro Dalcin +     pc - the preconditioner context
92916d9e3a6SLisandro Dalcin -     name - either  pilut, parasails, boomeramg, euclid
93016d9e3a6SLisandro Dalcin 
93116d9e3a6SLisandro Dalcin    Options Database Keys:
93216d9e3a6SLisandro Dalcin    -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
93316d9e3a6SLisandro Dalcin 
93416d9e3a6SLisandro Dalcin    Level: intermediate
93516d9e3a6SLisandro Dalcin 
93616d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
93716d9e3a6SLisandro Dalcin            PCHYPRE
93816d9e3a6SLisandro Dalcin 
93916d9e3a6SLisandro Dalcin @*/
9407087cfbeSBarry Smith PetscErrorCode  PCHYPRESetType(PC pc,const char name[])
94116d9e3a6SLisandro Dalcin {
9424ac538c5SBarry Smith   PetscErrorCode ierr;
94316d9e3a6SLisandro Dalcin 
94416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9450700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
94616d9e3a6SLisandro Dalcin   PetscValidCharPointer(name,2);
9474ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCHYPRESetType_C",(PC,const char[]),(pc,name));CHKERRQ(ierr);
94816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
94916d9e3a6SLisandro Dalcin }
95016d9e3a6SLisandro Dalcin 
95116d9e3a6SLisandro Dalcin #undef __FUNCT__
95216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType"
95316d9e3a6SLisandro Dalcin /*@C
95416d9e3a6SLisandro Dalcin      PCHYPREGetType - Gets which hypre preconditioner you are using
95516d9e3a6SLisandro Dalcin 
95616d9e3a6SLisandro Dalcin    Input Parameter:
95716d9e3a6SLisandro Dalcin .     pc - the preconditioner context
95816d9e3a6SLisandro Dalcin 
95916d9e3a6SLisandro Dalcin    Output Parameter:
96016d9e3a6SLisandro Dalcin .     name - either  pilut, parasails, boomeramg, euclid
96116d9e3a6SLisandro Dalcin 
96216d9e3a6SLisandro Dalcin    Level: intermediate
96316d9e3a6SLisandro Dalcin 
96416d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC,
96516d9e3a6SLisandro Dalcin            PCHYPRE
96616d9e3a6SLisandro Dalcin 
96716d9e3a6SLisandro Dalcin @*/
9687087cfbeSBarry Smith PetscErrorCode  PCHYPREGetType(PC pc,const char *name[])
96916d9e3a6SLisandro Dalcin {
9704ac538c5SBarry Smith   PetscErrorCode ierr;
97116d9e3a6SLisandro Dalcin 
97216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9730700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
97416d9e3a6SLisandro Dalcin   PetscValidPointer(name,2);
9754ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCHYPREGetType_C",(PC,const char*[]),(pc,name));CHKERRQ(ierr);
97616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
97716d9e3a6SLisandro Dalcin }
97816d9e3a6SLisandro Dalcin 
97916d9e3a6SLisandro Dalcin /*MC
98016d9e3a6SLisandro Dalcin      PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre
98116d9e3a6SLisandro Dalcin 
98216d9e3a6SLisandro Dalcin    Options Database Keys:
98316d9e3a6SLisandro Dalcin +   -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
98416d9e3a6SLisandro Dalcin -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
98516d9e3a6SLisandro Dalcin           preconditioner
98616d9e3a6SLisandro Dalcin 
98716d9e3a6SLisandro Dalcin    Level: intermediate
98816d9e3a6SLisandro Dalcin 
98916d9e3a6SLisandro Dalcin    Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()),
99016d9e3a6SLisandro Dalcin           the many hypre options can ONLY be set via the options database (e.g. the command line
99116d9e3a6SLisandro Dalcin           or with PetscOptionsSetValue(), there are no functions to set them)
99216d9e3a6SLisandro Dalcin 
99316d9e3a6SLisandro Dalcin           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations
9940f1074feSSatish Balay           (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if
9950f1074feSSatish Balay           -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner
9960f1074feSSatish Balay           (-pc_hypre_boomeramg_rtol should be set to 0.0 - the default - to strictly use a fixed number of
9978f87f92bSBarry Smith           iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations
9980f1074feSSatish Balay           and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10
9990f1074feSSatish Balay           then AT MOST twenty V-cycles of boomeramg will be called.
100016d9e3a6SLisandro Dalcin 
10010f1074feSSatish Balay            Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation
10020f1074feSSatish Balay            (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry.
10030f1074feSSatish Balay            Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi.
100416d9e3a6SLisandro Dalcin           If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly
100516d9e3a6SLisandro Dalcin           and use -ksp_max_it to control the number of V-cycles.
100616d9e3a6SLisandro Dalcin           (see the PETSc FAQ.html at the PETSc website under the Documentation tab).
100716d9e3a6SLisandro Dalcin 
100816d9e3a6SLisandro Dalcin           2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option
100916d9e3a6SLisandro Dalcin           -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L.
101016d9e3a6SLisandro Dalcin 
10119e5bc791SBarry Smith           See PCPFMG for access to the hypre Struct PFMG solver
10129e5bc791SBarry Smith 
101316d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
10149e5bc791SBarry Smith            PCHYPRESetType(), PCPFMG
101516d9e3a6SLisandro Dalcin 
101616d9e3a6SLisandro Dalcin M*/
101716d9e3a6SLisandro Dalcin 
101816d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
101916d9e3a6SLisandro Dalcin #undef __FUNCT__
102016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCCreate_HYPRE"
10217087cfbeSBarry Smith PetscErrorCode  PCCreate_HYPRE(PC pc)
102216d9e3a6SLisandro Dalcin {
102316d9e3a6SLisandro Dalcin   PC_HYPRE       *jac;
102416d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
102516d9e3a6SLisandro Dalcin 
102616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
102738f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_HYPRE,&jac);CHKERRQ(ierr);
102816d9e3a6SLisandro Dalcin   pc->data                 = jac;
102916d9e3a6SLisandro Dalcin   pc->ops->destroy         = PCDestroy_HYPRE;
103016d9e3a6SLisandro Dalcin   pc->ops->setfromoptions  = PCSetFromOptions_HYPRE;
103116d9e3a6SLisandro Dalcin   pc->ops->setup           = PCSetUp_HYPRE;
103216d9e3a6SLisandro Dalcin   pc->ops->apply           = PCApply_HYPRE;
103316d9e3a6SLisandro Dalcin   jac->comm_hypre          = MPI_COMM_NULL;
103416d9e3a6SLisandro Dalcin   jac->hypre_type          = PETSC_NULL;
103516d9e3a6SLisandro Dalcin   /* duplicate communicator for hypre */
10367adad957SLisandro Dalcin   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(jac->comm_hypre));CHKERRQ(ierr);
103716d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr);
103816d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr);
103916d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
104016d9e3a6SLisandro Dalcin }
104116d9e3a6SLisandro Dalcin EXTERN_C_END
1042ebc551c0SBarry Smith 
1043f91d8e95SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------*/
1044f91d8e95SBarry Smith 
1045b862ddfaSBarry Smith /* this include is needed ONLY to allow access to the private data inside the Mat object specific to hypre */
1046b45d2f2cSJed Brown #include <petsc-private/matimpl.h>
1047ebc551c0SBarry Smith 
1048ebc551c0SBarry Smith typedef struct {
104968326731SBarry Smith   MPI_Comm            hcomm;       /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */
1050f91d8e95SBarry Smith   HYPRE_StructSolver  hsolver;
10519e5bc791SBarry Smith 
10529e5bc791SBarry Smith   /* keep copy of PFMG options used so may view them */
10534ddd07fcSJed Brown   PetscInt            its;
10549e5bc791SBarry Smith   double              tol;
10554ddd07fcSJed Brown   PetscInt            relax_type;
10564ddd07fcSJed Brown   PetscInt            rap_type;
10574ddd07fcSJed Brown   PetscInt            num_pre_relax,num_post_relax;
10584ddd07fcSJed Brown   PetscInt            max_levels;
1059ebc551c0SBarry Smith } PC_PFMG;
1060ebc551c0SBarry Smith 
1061ebc551c0SBarry Smith #undef __FUNCT__
1062ebc551c0SBarry Smith #define __FUNCT__ "PCDestroy_PFMG"
1063ebc551c0SBarry Smith PetscErrorCode PCDestroy_PFMG(PC pc)
1064ebc551c0SBarry Smith {
1065ebc551c0SBarry Smith   PetscErrorCode ierr;
1066f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1067ebc551c0SBarry Smith 
1068ebc551c0SBarry Smith   PetscFunctionBegin;
1069fd3f9acdSBarry Smith   if (ex->hsolver) {PetscStackCallStandard(HYPRE_StructPFMGDestroy,(ex->hsolver));}
1070f91d8e95SBarry Smith   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1071c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1072ebc551c0SBarry Smith   PetscFunctionReturn(0);
1073ebc551c0SBarry Smith }
1074ebc551c0SBarry Smith 
10759e5bc791SBarry Smith static const char *PFMGRelaxType[]   = {"Jacobi","Weighted-Jacobi","symmetric-Red/Black-Gauss-Seidel","Red/Black-Gauss-Seidel"};
10769e5bc791SBarry Smith static const char *PFMGRAPType[]   = {"Galerkin","non-Galerkin"};
10779e5bc791SBarry Smith 
1078ebc551c0SBarry Smith #undef __FUNCT__
1079ebc551c0SBarry Smith #define __FUNCT__ "PCView_PFMG"
1080ebc551c0SBarry Smith PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer)
1081ebc551c0SBarry Smith {
1082ebc551c0SBarry Smith   PetscErrorCode ierr;
1083ace3abfcSBarry Smith   PetscBool      iascii;
1084f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1085ebc551c0SBarry Smith 
1086ebc551c0SBarry Smith   PetscFunctionBegin;
1087251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10889e5bc791SBarry Smith   if (iascii) {
10899e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG preconditioning\n");CHKERRQ(ierr);
10909e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
10919e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
10929e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
10939e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: RAP type %s\n",PFMGRAPType[ex->rap_type]);CHKERRQ(ierr);
10949e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
10953b46a515SGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max levels %d\n",ex->max_levels);CHKERRQ(ierr);
10969e5bc791SBarry Smith   }
1097ebc551c0SBarry Smith   PetscFunctionReturn(0);
1098ebc551c0SBarry Smith }
1099ebc551c0SBarry Smith 
11009e5bc791SBarry Smith 
1101ebc551c0SBarry Smith #undef __FUNCT__
1102ebc551c0SBarry Smith #define __FUNCT__ "PCSetFromOptions_PFMG"
1103ebc551c0SBarry Smith PetscErrorCode PCSetFromOptions_PFMG(PC pc)
1104ebc551c0SBarry Smith {
1105ebc551c0SBarry Smith   PetscErrorCode ierr;
1106f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1107ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
1108ebc551c0SBarry Smith 
1109ebc551c0SBarry Smith   PetscFunctionBegin;
1110ebc551c0SBarry Smith   ierr = PetscOptionsHead("PFMG options");CHKERRQ(ierr);
1111acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
111268326731SBarry Smith   if (flg) {
1113a0324ebeSBarry Smith     int level=3;
1114fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_StructPFMGSetPrintLevel,(ex->hsolver,level));
111568326731SBarry Smith   }
11169e5bc791SBarry Smith   ierr = PetscOptionsInt("-pc_pfmg_its","Number of iterations of PFMG to use as preconditioner","HYPRE_StructPFMGSetMaxIter",ex->its,&ex->its,PETSC_NULL);CHKERRQ(ierr);
1117fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(ex->hsolver,ex->its));
11189e5bc791SBarry 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,PETSC_NULL);CHKERRQ(ierr);
1119fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetNumPreRelax,(ex->hsolver,ex->num_pre_relax));
11209e5bc791SBarry 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,PETSC_NULL);CHKERRQ(ierr);
1121fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetNumPostRelax,(ex->hsolver,ex->num_post_relax));
11229e5bc791SBarry Smith 
11233b46a515SGlenn Hammond   ierr = PetscOptionsInt("-pc_pfmg_max_levels","Max Levels for MG hierarchy","HYPRE_StructPFMGSetMaxLevels",ex->max_levels,&ex->max_levels,PETSC_NULL);CHKERRQ(ierr);
1124fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetMaxLevels,(ex->hsolver,ex->max_levels));
11253b46a515SGlenn Hammond 
11269e5bc791SBarry Smith   ierr = PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr);
1127fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(ex->hsolver,ex->tol));
1128a669f990SJed Brown   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,PETSC_NULL);CHKERRQ(ierr);
1129fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetRelaxType,(ex->hsolver, ex->relax_type));
1130a669f990SJed Brown   ierr = PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,ALEN(PFMGRAPType),PFMGRAPType[ex->rap_type],&ex->rap_type,PETSC_NULL);CHKERRQ(ierr);
1131fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetRAPType,(ex->hsolver, ex->rap_type));
1132ebc551c0SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1133ebc551c0SBarry Smith   PetscFunctionReturn(0);
1134ebc551c0SBarry Smith }
1135ebc551c0SBarry Smith 
1136f91d8e95SBarry Smith #undef __FUNCT__
1137f91d8e95SBarry Smith #define __FUNCT__ "PCApply_PFMG"
1138f91d8e95SBarry Smith PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y)
1139f91d8e95SBarry Smith {
1140f91d8e95SBarry Smith   PetscErrorCode  ierr;
1141f91d8e95SBarry Smith   PC_PFMG         *ex = (PC_PFMG*) pc->data;
1142f91d8e95SBarry Smith   PetscScalar     *xx,*yy;
11434ddd07fcSJed Brown   PetscInt        ilower[3],iupper[3];
114468326731SBarry Smith   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
1145f91d8e95SBarry Smith 
1146f91d8e95SBarry Smith   PetscFunctionBegin;
1147aa219208SBarry Smith   ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1148f91d8e95SBarry Smith   iupper[0] += ilower[0] - 1;
1149f91d8e95SBarry Smith   iupper[1] += ilower[1] - 1;
1150f91d8e95SBarry Smith   iupper[2] += ilower[2] - 1;
1151f91d8e95SBarry Smith 
1152f91d8e95SBarry Smith   /* copy x values over to hypre */
1153fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorSetConstantValues,(mx->hb,0.0));
1154f91d8e95SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1155fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorSetBoxValues,(mx->hb,ilower,iupper,xx));
1156f91d8e95SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1157fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorAssemble,(mx->hb));
1158fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSolve,(ex->hsolver,mx->hmat,mx->hb,mx->hx));
1159f91d8e95SBarry Smith 
1160f91d8e95SBarry Smith   /* copy solution values back to PETSc */
1161f91d8e95SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1162fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorGetBoxValues,(mx->hx,ilower,iupper,yy));
1163f91d8e95SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1164f91d8e95SBarry Smith   PetscFunctionReturn(0);
1165f91d8e95SBarry Smith }
1166f91d8e95SBarry Smith 
11679e5bc791SBarry Smith #undef __FUNCT__
11689e5bc791SBarry Smith #define __FUNCT__ "PCApplyRichardson_PFMG"
1169ace3abfcSBarry 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)
11709e5bc791SBarry Smith {
11719e5bc791SBarry Smith   PC_PFMG        *jac = (PC_PFMG*)pc->data;
11729e5bc791SBarry Smith   PetscErrorCode ierr;
11734ddd07fcSJed Brown   PetscInt       oits;
11749e5bc791SBarry Smith 
11759e5bc791SBarry Smith   PetscFunctionBegin;
1176fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(jac->hsolver,its*jac->its));
1177fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(jac->hsolver,rtol));
11789e5bc791SBarry Smith 
11799e5bc791SBarry Smith   ierr = PCApply_PFMG(pc,b,y);CHKERRQ(ierr);
1180fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGGetNumIterations,(jac->hsolver,&oits));
11819e5bc791SBarry Smith   *outits = oits;
11829e5bc791SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
11839e5bc791SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
1184fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(jac->hsolver,jac->tol));
1185fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(jac->hsolver,jac->its));
11869e5bc791SBarry Smith   PetscFunctionReturn(0);
11879e5bc791SBarry Smith }
11889e5bc791SBarry Smith 
11899e5bc791SBarry Smith 
11903a32d3dbSGlenn Hammond #undef __FUNCT__
11913a32d3dbSGlenn Hammond #define __FUNCT__ "PCSetUp_PFMG"
11923a32d3dbSGlenn Hammond PetscErrorCode PCSetUp_PFMG(PC pc)
11933a32d3dbSGlenn Hammond {
11943a32d3dbSGlenn Hammond   PetscErrorCode  ierr;
11953a32d3dbSGlenn Hammond   PC_PFMG         *ex = (PC_PFMG*) pc->data;
11963a32d3dbSGlenn Hammond   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
1197ace3abfcSBarry Smith   PetscBool       flg;
11983a32d3dbSGlenn Hammond 
11993a32d3dbSGlenn Hammond   PetscFunctionBegin;
1200251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESTRUCT,&flg);CHKERRQ(ierr);
1201e7e72b3dSBarry Smith   if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESTRUCT with this preconditioner");
12023a32d3dbSGlenn Hammond 
12033a32d3dbSGlenn Hammond   /* create the hypre solver object and set its information */
12043a32d3dbSGlenn Hammond   if (ex->hsolver) {
1205fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_StructPFMGDestroy,(ex->hsolver));
12063a32d3dbSGlenn Hammond   }
1207fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
12083a32d3dbSGlenn Hammond   ierr = PCSetFromOptions_PFMG(pc);CHKERRQ(ierr);
1209fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetup,(ex->hsolver,mx->hmat,mx->hb,mx->hx));
1210fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetZeroGuess,(ex->hsolver));
12113a32d3dbSGlenn Hammond   PetscFunctionReturn(0);
12123a32d3dbSGlenn Hammond }
12133a32d3dbSGlenn Hammond 
1214ebc551c0SBarry Smith 
1215ebc551c0SBarry Smith /*MC
1216ebc551c0SBarry Smith      PCPFMG - the hypre PFMG multigrid solver
1217ebc551c0SBarry Smith 
1218ebc551c0SBarry Smith    Level: advanced
1219ebc551c0SBarry Smith 
12209e5bc791SBarry Smith    Options Database:
12219e5bc791SBarry Smith + -pc_pfmg_its <its> number of iterations of PFMG to use as preconditioner
12229e5bc791SBarry Smith . -pc_pfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
12239e5bc791SBarry Smith . -pc_pfmg_num_post_relax <steps> number of smoothing steps after coarse grid
12249e5bc791SBarry Smith . -pc_pfmg_tol <tol> tolerance of PFMG
12259e5bc791SBarry 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
12269e5bc791SBarry Smith - -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin
1227f91d8e95SBarry Smith 
12289e5bc791SBarry Smith    Notes:  This is for CELL-centered descretizations
12299e5bc791SBarry Smith 
12308e395302SJed Brown            This must be used with the MATHYPRESTRUCT matrix type.
1231aa219208SBarry Smith            This is less general than in hypre, it supports only one block per process defined by a PETSc DMDA.
12329e5bc791SBarry Smith 
12339e5bc791SBarry Smith .seealso:  PCMG, MATHYPRESTRUCT
1234ebc551c0SBarry Smith M*/
1235ebc551c0SBarry Smith 
1236ebc551c0SBarry Smith EXTERN_C_BEGIN
1237ebc551c0SBarry Smith #undef __FUNCT__
1238ebc551c0SBarry Smith #define __FUNCT__ "PCCreate_PFMG"
12397087cfbeSBarry Smith PetscErrorCode  PCCreate_PFMG(PC pc)
1240ebc551c0SBarry Smith {
1241ebc551c0SBarry Smith   PetscErrorCode ierr;
1242ebc551c0SBarry Smith   PC_PFMG        *ex;
1243ebc551c0SBarry Smith 
1244ebc551c0SBarry Smith   PetscFunctionBegin;
1245ebc551c0SBarry Smith   ierr = PetscNew(PC_PFMG,&ex);CHKERRQ(ierr);\
124668326731SBarry Smith   pc->data = ex;
1247ebc551c0SBarry Smith 
12489e5bc791SBarry Smith   ex->its            = 1;
12499e5bc791SBarry Smith   ex->tol            = 1.e-8;
12509e5bc791SBarry Smith   ex->relax_type     = 1;
12519e5bc791SBarry Smith   ex->rap_type       = 0;
12529e5bc791SBarry Smith   ex->num_pre_relax  = 1;
12539e5bc791SBarry Smith   ex->num_post_relax = 1;
12543b46a515SGlenn Hammond   ex->max_levels     = 0;
12559e5bc791SBarry Smith 
1256ebc551c0SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_PFMG;
1257ebc551c0SBarry Smith   pc->ops->view            = PCView_PFMG;
1258ebc551c0SBarry Smith   pc->ops->destroy         = PCDestroy_PFMG;
1259f91d8e95SBarry Smith   pc->ops->apply           = PCApply_PFMG;
12609e5bc791SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_PFMG;
126168326731SBarry Smith   pc->ops->setup           = PCSetUp_PFMG;
1262f91d8e95SBarry Smith   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr);
1263fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
1264ebc551c0SBarry Smith   PetscFunctionReturn(0);
1265ebc551c0SBarry Smith }
1266ebc551c0SBarry Smith EXTERN_C_END
1267d851a50bSGlenn Hammond 
1268325fc9f4SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------------------------*/
1269325fc9f4SBarry Smith 
1270d851a50bSGlenn Hammond /* we know we are working with a HYPRE_SStructMatrix */
1271d851a50bSGlenn Hammond typedef struct {
1272d851a50bSGlenn Hammond   MPI_Comm            hcomm;       /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */
1273d851a50bSGlenn Hammond   HYPRE_SStructSolver  ss_solver;
1274d851a50bSGlenn Hammond 
1275d851a50bSGlenn Hammond   /* keep copy of SYSPFMG options used so may view them */
12764ddd07fcSJed Brown   PetscInt            its;
1277d851a50bSGlenn Hammond   double              tol;
12784ddd07fcSJed Brown   PetscInt            relax_type;
12794ddd07fcSJed Brown   PetscInt            num_pre_relax,num_post_relax;
1280d851a50bSGlenn Hammond } PC_SysPFMG;
1281d851a50bSGlenn Hammond 
1282d851a50bSGlenn Hammond #undef __FUNCT__
1283d851a50bSGlenn Hammond #define __FUNCT__ "PCDestroy_SysPFMG"
1284d851a50bSGlenn Hammond PetscErrorCode PCDestroy_SysPFMG(PC pc)
1285d851a50bSGlenn Hammond {
1286d851a50bSGlenn Hammond   PetscErrorCode ierr;
1287d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1288d851a50bSGlenn Hammond 
1289d851a50bSGlenn Hammond   PetscFunctionBegin;
1290fd3f9acdSBarry Smith   if (ex->ss_solver) {PetscStackCallStandard(HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));}
1291d851a50bSGlenn Hammond   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1292c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1293d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1294d851a50bSGlenn Hammond }
1295d851a50bSGlenn Hammond 
1296d851a50bSGlenn Hammond static const char *SysPFMGRelaxType[]   = {"Weighted-Jacobi","Red/Black-Gauss-Seidel"};
1297d851a50bSGlenn Hammond 
1298d851a50bSGlenn Hammond #undef __FUNCT__
1299d851a50bSGlenn Hammond #define __FUNCT__ "PCView_SysPFMG"
1300d851a50bSGlenn Hammond PetscErrorCode PCView_SysPFMG(PC pc,PetscViewer viewer)
1301d851a50bSGlenn Hammond {
1302d851a50bSGlenn Hammond   PetscErrorCode ierr;
1303ace3abfcSBarry Smith   PetscBool      iascii;
1304d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1305d851a50bSGlenn Hammond 
1306d851a50bSGlenn Hammond   PetscFunctionBegin;
1307251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1308d851a50bSGlenn Hammond   if (iascii) {
1309d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG preconditioning\n");CHKERRQ(ierr);
1310d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
1311d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
1312d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
1313d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
1314d851a50bSGlenn Hammond   }
1315d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1316d851a50bSGlenn Hammond }
1317d851a50bSGlenn Hammond 
1318d851a50bSGlenn Hammond 
1319d851a50bSGlenn Hammond #undef __FUNCT__
1320d851a50bSGlenn Hammond #define __FUNCT__ "PCSetFromOptions_SysPFMG"
1321d851a50bSGlenn Hammond PetscErrorCode PCSetFromOptions_SysPFMG(PC pc)
1322d851a50bSGlenn Hammond {
1323d851a50bSGlenn Hammond   PetscErrorCode ierr;
1324d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1325ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
1326d851a50bSGlenn Hammond 
1327d851a50bSGlenn Hammond   PetscFunctionBegin;
1328d851a50bSGlenn Hammond   ierr = PetscOptionsHead("SysPFMG options");CHKERRQ(ierr);
1329acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_syspfmg_print_statistics","Print statistics","HYPRE_SStructSysPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1330d851a50bSGlenn Hammond   if (flg) {
1331d851a50bSGlenn Hammond     int level=3;
1332fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructSysPFMGSetPrintLevel,(ex->ss_solver,level));
1333d851a50bSGlenn Hammond   }
1334d851a50bSGlenn Hammond   ierr = PetscOptionsInt("-pc_syspfmg_its","Number of iterations of SysPFMG to use as preconditioner","HYPRE_SStructSysPFMGSetMaxIter",ex->its,&ex->its,PETSC_NULL);CHKERRQ(ierr);
1335fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(ex->ss_solver,ex->its));
1336d851a50bSGlenn Hammond   ierr = PetscOptionsInt("-pc_syspfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_SStructSysPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,PETSC_NULL);CHKERRQ(ierr);
1337fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetNumPreRelax,(ex->ss_solver,ex->num_pre_relax));
1338d851a50bSGlenn Hammond   ierr = PetscOptionsInt("-pc_syspfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_SStructSysPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,PETSC_NULL);CHKERRQ(ierr);
1339fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetNumPostRelax,(ex->ss_solver,ex->num_post_relax));
1340d851a50bSGlenn Hammond 
1341d851a50bSGlenn Hammond   ierr = PetscOptionsReal("-pc_syspfmg_tol","Tolerance of SysPFMG","HYPRE_SStructSysPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr);
1342fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(ex->ss_solver,ex->tol));
1343d851a50bSGlenn Hammond   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,PETSC_NULL);CHKERRQ(ierr);
1344fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetRelaxType,(ex->ss_solver, ex->relax_type));
1345d851a50bSGlenn Hammond   ierr = PetscOptionsTail();CHKERRQ(ierr);
1346d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1347d851a50bSGlenn Hammond }
1348d851a50bSGlenn Hammond 
1349d851a50bSGlenn Hammond #undef __FUNCT__
1350d851a50bSGlenn Hammond #define __FUNCT__ "PCApply_SysPFMG"
1351d851a50bSGlenn Hammond PetscErrorCode PCApply_SysPFMG(PC pc,Vec x,Vec y)
1352d851a50bSGlenn Hammond {
1353d851a50bSGlenn Hammond   PetscErrorCode    ierr;
1354d851a50bSGlenn Hammond   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
1355d851a50bSGlenn Hammond   PetscScalar      *xx,*yy;
13564ddd07fcSJed Brown   PetscInt          ilower[3],iupper[3];
1357d851a50bSGlenn Hammond   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct *)(pc->pmat->data);
13584ddd07fcSJed Brown   PetscInt          ordering= mx->dofs_order;
13594ddd07fcSJed Brown   PetscInt          nvars= mx->nvars;
13604ddd07fcSJed Brown   PetscInt          part= 0;
13614ddd07fcSJed Brown   PetscInt          size;
13624ddd07fcSJed Brown   PetscInt          i;
1363d851a50bSGlenn Hammond 
1364d851a50bSGlenn Hammond   PetscFunctionBegin;
1365aa219208SBarry Smith   ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1366d851a50bSGlenn Hammond   iupper[0] += ilower[0] - 1;
1367d851a50bSGlenn Hammond   iupper[1] += ilower[1] - 1;
1368d851a50bSGlenn Hammond   iupper[2] += ilower[2] - 1;
1369d851a50bSGlenn Hammond 
1370d851a50bSGlenn Hammond   size= 1;
1371d851a50bSGlenn Hammond   for (i= 0; i< 3; i++) {
1372d851a50bSGlenn Hammond      size*= (iupper[i]-ilower[i]+1);
1373d851a50bSGlenn Hammond   }
1374d851a50bSGlenn Hammond   /* copy x values over to hypre for variable ordering */
1375d851a50bSGlenn Hammond   if (ordering) {
1376fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
1377d851a50bSGlenn Hammond      ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1378d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
1379fd3f9acdSBarry Smith        PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,ilower,iupper,i,xx+(size*i)));
1380d851a50bSGlenn Hammond      }
1381d851a50bSGlenn Hammond      ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1382fd3f9acdSBarry Smith      PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
1383fd3f9acdSBarry Smith      PetscStackCallStandard(HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
1384fd3f9acdSBarry Smith      PetscStackCallStandard(HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1385d851a50bSGlenn Hammond 
1386d851a50bSGlenn Hammond      /* copy solution values back to PETSc */
1387d851a50bSGlenn Hammond      ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1388d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
1389fd3f9acdSBarry Smith        PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,ilower,iupper,i,yy+(size*i)));
1390d851a50bSGlenn Hammond      }
1391d851a50bSGlenn Hammond      ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1392a65764d7SBarry Smith   } else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
1393d851a50bSGlenn Hammond      PetscScalar     *z;
13944ddd07fcSJed Brown      PetscInt         j, k;
1395d851a50bSGlenn Hammond 
1396d851a50bSGlenn Hammond      ierr = PetscMalloc(nvars*size*sizeof(PetscScalar),&z);CHKERRQ(ierr);
1397fd3f9acdSBarry Smith      PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
1398d851a50bSGlenn Hammond      ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1399d851a50bSGlenn Hammond 
1400d851a50bSGlenn Hammond      /* transform nodal to hypre's variable ordering for sys_pfmg */
1401d851a50bSGlenn Hammond      for (i= 0; i< size; i++) {
1402d851a50bSGlenn Hammond         k= i*nvars;
1403d851a50bSGlenn Hammond         for (j= 0; j< nvars; j++) {
1404d851a50bSGlenn Hammond            z[j*size+i]= xx[k+j];
1405d851a50bSGlenn Hammond         }
1406d851a50bSGlenn Hammond      }
1407d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
1408fd3f9acdSBarry Smith        PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,ilower,iupper,i,z+(size*i)));
1409d851a50bSGlenn Hammond      }
1410d851a50bSGlenn Hammond      ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1411fd3f9acdSBarry Smith      PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
1412fd3f9acdSBarry Smith      PetscStackCallStandard(HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1413d851a50bSGlenn Hammond 
1414d851a50bSGlenn Hammond      /* copy solution values back to PETSc */
1415d851a50bSGlenn Hammond      ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1416d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
1417fd3f9acdSBarry Smith        PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,ilower,iupper,i,z+(size*i)));
1418d851a50bSGlenn Hammond      }
1419d851a50bSGlenn Hammond      /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
1420d851a50bSGlenn Hammond      for (i= 0; i< size; i++) {
1421d851a50bSGlenn Hammond         k= i*nvars;
1422d851a50bSGlenn Hammond         for (j= 0; j< nvars; j++) {
1423d851a50bSGlenn Hammond            yy[k+j]= z[j*size+i];
1424d851a50bSGlenn Hammond         }
1425d851a50bSGlenn Hammond      }
1426d851a50bSGlenn Hammond      ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1427d851a50bSGlenn Hammond      ierr = PetscFree(z);CHKERRQ(ierr);
1428d851a50bSGlenn Hammond   }
1429d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1430d851a50bSGlenn Hammond }
1431d851a50bSGlenn Hammond 
1432d851a50bSGlenn Hammond #undef __FUNCT__
1433d851a50bSGlenn Hammond #define __FUNCT__ "PCApplyRichardson_SysPFMG"
1434ace3abfcSBarry 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)
1435d851a50bSGlenn Hammond {
1436d851a50bSGlenn Hammond   PC_SysPFMG    *jac = (PC_SysPFMG*)pc->data;
1437d851a50bSGlenn Hammond   PetscErrorCode ierr;
14384ddd07fcSJed Brown   PetscInt       oits;
1439d851a50bSGlenn Hammond 
1440d851a50bSGlenn Hammond   PetscFunctionBegin;
1441fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,its*jac->its));
1442fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,rtol));
1443d851a50bSGlenn Hammond   ierr = PCApply_SysPFMG(pc,b,y);CHKERRQ(ierr);
1444fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGGetNumIterations,(jac->ss_solver,&oits));
1445d851a50bSGlenn Hammond   *outits = oits;
1446d851a50bSGlenn Hammond   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
1447d851a50bSGlenn Hammond   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
1448fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,jac->tol));
1449fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,jac->its));
1450d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1451d851a50bSGlenn Hammond }
1452d851a50bSGlenn Hammond 
1453d851a50bSGlenn Hammond 
1454d851a50bSGlenn Hammond #undef __FUNCT__
1455d851a50bSGlenn Hammond #define __FUNCT__ "PCSetUp_SysPFMG"
1456d851a50bSGlenn Hammond PetscErrorCode PCSetUp_SysPFMG(PC pc)
1457d851a50bSGlenn Hammond {
1458d851a50bSGlenn Hammond   PetscErrorCode    ierr;
1459d851a50bSGlenn Hammond   PC_SysPFMG        *ex = (PC_SysPFMG*) pc->data;
1460d851a50bSGlenn Hammond   Mat_HYPRESStruct  *mx = (Mat_HYPRESStruct *)(pc->pmat->data);
1461ace3abfcSBarry Smith   PetscBool         flg;
1462d851a50bSGlenn Hammond 
1463d851a50bSGlenn Hammond   PetscFunctionBegin;
1464251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESSTRUCT,&flg);CHKERRQ(ierr);
1465e7e72b3dSBarry Smith   if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESSTRUCT with this preconditioner");
1466d851a50bSGlenn Hammond 
1467d851a50bSGlenn Hammond   /* create the hypre sstruct solver object and set its information */
1468d851a50bSGlenn Hammond   if (ex->ss_solver) {
1469fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));
1470d851a50bSGlenn Hammond   }
1471fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
1472d851a50bSGlenn Hammond   ierr = PCSetFromOptions_SysPFMG(pc);CHKERRQ(ierr);
1473fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetZeroGuess,(ex->ss_solver));
1474fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetup,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1475d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1476d851a50bSGlenn Hammond }
1477d851a50bSGlenn Hammond 
1478d851a50bSGlenn Hammond 
1479d851a50bSGlenn Hammond /*MC
1480d851a50bSGlenn Hammond      PCSysPFMG - the hypre SysPFMG multigrid solver
1481d851a50bSGlenn Hammond 
1482d851a50bSGlenn Hammond    Level: advanced
1483d851a50bSGlenn Hammond 
1484d851a50bSGlenn Hammond    Options Database:
1485d851a50bSGlenn Hammond + -pc_syspfmg_its <its> number of iterations of SysPFMG to use as preconditioner
1486d851a50bSGlenn Hammond . -pc_syspfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
1487d851a50bSGlenn Hammond . -pc_syspfmg_num_post_relax <steps> number of smoothing steps after coarse grid
1488d851a50bSGlenn Hammond . -pc_syspfmg_tol <tol> tolerance of SysPFMG
1489d851a50bSGlenn Hammond . -pc_syspfmg_relax_type -relaxation type for the up and down cycles, one of Weighted-Jacobi,Red/Black-Gauss-Seidel
1490d851a50bSGlenn Hammond 
1491d851a50bSGlenn Hammond    Notes:  This is for CELL-centered descretizations
1492d851a50bSGlenn Hammond 
1493f6680f47SSatish Balay            This must be used with the MATHYPRESSTRUCT matrix type.
1494aa219208SBarry Smith            This is less general than in hypre, it supports only one part, and one block per process defined by a PETSc DMDA.
1495d851a50bSGlenn Hammond            Also, only cell-centered variables.
1496d851a50bSGlenn Hammond 
1497d851a50bSGlenn Hammond .seealso:  PCMG, MATHYPRESSTRUCT
1498d851a50bSGlenn Hammond M*/
1499d851a50bSGlenn Hammond 
1500d851a50bSGlenn Hammond EXTERN_C_BEGIN
1501d851a50bSGlenn Hammond #undef __FUNCT__
1502d851a50bSGlenn Hammond #define __FUNCT__ "PCCreate_SysPFMG"
15037087cfbeSBarry Smith PetscErrorCode  PCCreate_SysPFMG(PC pc)
1504d851a50bSGlenn Hammond {
1505d851a50bSGlenn Hammond   PetscErrorCode     ierr;
1506d851a50bSGlenn Hammond   PC_SysPFMG        *ex;
1507d851a50bSGlenn Hammond 
1508d851a50bSGlenn Hammond   PetscFunctionBegin;
1509d851a50bSGlenn Hammond   ierr = PetscNew(PC_SysPFMG,&ex);CHKERRQ(ierr);\
1510d851a50bSGlenn Hammond   pc->data = ex;
1511d851a50bSGlenn Hammond 
1512d851a50bSGlenn Hammond   ex->its            = 1;
1513d851a50bSGlenn Hammond   ex->tol            = 1.e-8;
1514d851a50bSGlenn Hammond   ex->relax_type     = 1;
1515d851a50bSGlenn Hammond   ex->num_pre_relax  = 1;
1516d851a50bSGlenn Hammond   ex->num_post_relax = 1;
1517d851a50bSGlenn Hammond 
1518d851a50bSGlenn Hammond   pc->ops->setfromoptions  = PCSetFromOptions_SysPFMG;
1519d851a50bSGlenn Hammond   pc->ops->view            = PCView_SysPFMG;
1520d851a50bSGlenn Hammond   pc->ops->destroy         = PCDestroy_SysPFMG;
1521d851a50bSGlenn Hammond   pc->ops->apply           = PCApply_SysPFMG;
1522d851a50bSGlenn Hammond   pc->ops->applyrichardson = PCApplyRichardson_SysPFMG;
1523d851a50bSGlenn Hammond   pc->ops->setup           = PCSetUp_SysPFMG;
1524d851a50bSGlenn Hammond   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr);
1525fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
1526d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1527d851a50bSGlenn Hammond }
1528d851a50bSGlenn Hammond EXTERN_C_END
1529