xref: /petsc/src/ksp/pc/impls/hypre/hypre.c (revision 65de4495e32eb7796777f05a3f700d2e59d93713)
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. */
10530e6f737SJed Brown     PetscStackCallHypre(0,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) {
12530e6f737SJed Brown       PetscStackCallHypre(0,HYPRE_BoomerAMGSetNumFunctions,(jac->hsolver,bs));
12616d9e3a6SLisandro Dalcin     }
12716d9e3a6SLisandro Dalcin   };
12816d9e3a6SLisandro Dalcin   ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr);
12930e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
13030e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJVectorGetObject,(jac->b,(void**)&bv));
13130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJVectorGetObject,(jac->x,(void**)&xv));
13230e6f737SJed Brown   PetscStackCallHypre("HYPRE_SetupXXX",(*jac->setup),(jac->hsolver,hmat,bv,xv));
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;\
14516d9e3a6SLisandro Dalcin    local_vector->data = newvalue;}
14616d9e3a6SLisandro Dalcin 
14716d9e3a6SLisandro Dalcin #undef __FUNCT__
14816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApply_HYPRE"
14916d9e3a6SLisandro Dalcin static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
15016d9e3a6SLisandro Dalcin {
15116d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
15216d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
15316d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
15416d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
15516d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
15616d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
1574ddd07fcSJed Brown   PetscInt           hierr;
15816d9e3a6SLisandro Dalcin 
15916d9e3a6SLisandro Dalcin   PetscFunctionBegin;
16016d9e3a6SLisandro Dalcin   if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);}
16116d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
16216d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
16316d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
16416d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
16516d9e3a6SLisandro Dalcin 
16630e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
16730e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJVectorGetObject,(jac->b,(void**)&jbv));
16830e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJVectorGetObject,(jac->x,(void**)&jxv));
16916d9e3a6SLisandro Dalcin   hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv);
17065e19b50SBarry Smith   if (hierr && hierr != HYPRE_ERROR_CONV) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
17116d9e3a6SLisandro Dalcin   if (hierr) hypre__global_error = 0;
17216d9e3a6SLisandro Dalcin 
17316d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
17416d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
17516d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
17616d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
17716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
17816d9e3a6SLisandro Dalcin }
17916d9e3a6SLisandro Dalcin 
18016d9e3a6SLisandro Dalcin #undef __FUNCT__
18116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCDestroy_HYPRE"
18216d9e3a6SLisandro Dalcin static PetscErrorCode PCDestroy_HYPRE(PC pc)
18316d9e3a6SLisandro Dalcin {
18416d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
18516d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
18616d9e3a6SLisandro Dalcin 
18716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
18830e6f737SJed Brown   if (jac->ij) PetscStackCallHypre(0,HYPRE_IJMatrixDestroy,(jac->ij));
18930e6f737SJed Brown   if (jac->b)  PetscStackCallHypre(0,HYPRE_IJVectorDestroy,(jac->b));
19030e6f737SJed Brown   if (jac->x)  PetscStackCallHypre(0,HYPRE_IJVectorDestroy,(jac->x));
19130e6f737SJed Brown   if (jac->destroy) PetscStackCallHypre("HYPRE_DistroyXXX",(*jac->destroy),(jac->hsolver));
192503cfb0cSBarry Smith   ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr);
19316d9e3a6SLisandro Dalcin   if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);}
194c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
19516d9e3a6SLisandro Dalcin 
19616d9e3a6SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr);
19716d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);CHKERRQ(ierr);
19816d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);CHKERRQ(ierr);
19916d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
20016d9e3a6SLisandro Dalcin }
20116d9e3a6SLisandro Dalcin 
20216d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
20316d9e3a6SLisandro Dalcin #undef __FUNCT__
20416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut"
20516d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc)
20616d9e3a6SLisandro Dalcin {
20716d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
20816d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
209ace3abfcSBarry Smith   PetscBool      flag;
21016d9e3a6SLisandro Dalcin 
21116d9e3a6SLisandro Dalcin   PetscFunctionBegin;
21216d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr);
21316d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr);
21430e6f737SJed Brown   if (flag) PetscStackCallHypre(0,HYPRE_ParCSRPilutSetMaxIter,(jac->hsolver,jac->maxiter));
21516d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr);
21630e6f737SJed Brown   if (flag) PetscStackCallHypre(0,HYPRE_ParCSRPilutSetDropTolerance,(jac->hsolver,jac->tol));
21716d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr);
21830e6f737SJed Brown   if (flag) PetscStackCallHypre(0,HYPRE_ParCSRPilutSetFactorRowSize,(jac->hsolver,jac->factorrowsize));
21916d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
22016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
22116d9e3a6SLisandro Dalcin }
22216d9e3a6SLisandro Dalcin 
22316d9e3a6SLisandro Dalcin #undef __FUNCT__
22416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Pilut"
22516d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
22616d9e3a6SLisandro Dalcin {
22716d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
22816d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
229ace3abfcSBarry Smith   PetscBool      iascii;
23016d9e3a6SLisandro Dalcin 
23116d9e3a6SLisandro Dalcin   PetscFunctionBegin;
232251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
23316d9e3a6SLisandro Dalcin   if (iascii) {
23416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");CHKERRQ(ierr);
23516d9e3a6SLisandro Dalcin     if (jac->maxiter != PETSC_DEFAULT) {
23616d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr);
23716d9e3a6SLisandro Dalcin     } else {
23816d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr);
23916d9e3a6SLisandro Dalcin     }
24016d9e3a6SLisandro Dalcin     if (jac->tol != PETSC_DEFAULT) {
24116d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr);
24216d9e3a6SLisandro Dalcin     } else {
24316d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr);
24416d9e3a6SLisandro Dalcin     }
24516d9e3a6SLisandro Dalcin     if (jac->factorrowsize != PETSC_DEFAULT) {
24616d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr);
24716d9e3a6SLisandro Dalcin     } else {
24816d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");CHKERRQ(ierr);
24916d9e3a6SLisandro Dalcin     }
25016d9e3a6SLisandro Dalcin   }
25116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
25216d9e3a6SLisandro Dalcin }
25316d9e3a6SLisandro Dalcin 
25416d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
25516d9e3a6SLisandro Dalcin #undef __FUNCT__
25616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid"
25716d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc)
25816d9e3a6SLisandro Dalcin {
25916d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
26016d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
261ace3abfcSBarry Smith   PetscBool      flag;
262390e7148SBarry Smith   char           *args[8],levels[16];
263390e7148SBarry Smith   PetscInt       cnt = 0;
26416d9e3a6SLisandro Dalcin 
26516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
26616d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr);
26716d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr);
26816d9e3a6SLisandro Dalcin   if (flag) {
26965e19b50SBarry Smith     if (jac->levels < 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels);
2704ddd07fcSJed Brown     ierr = PetscSNPrintf(levels,sizeof levels,"%D",jac->levels);CHKERRQ(ierr);
271390e7148SBarry Smith     args[cnt++] = (char*)"-level"; args[cnt++] = levels;
27216d9e3a6SLisandro Dalcin   }
273acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);CHKERRQ(ierr);
27416d9e3a6SLisandro Dalcin   if (jac->bjilu) {
275390e7148SBarry Smith     args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1";
27616d9e3a6SLisandro Dalcin   }
27716d9e3a6SLisandro Dalcin 
278acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);CHKERRQ(ierr);
27916d9e3a6SLisandro Dalcin   if (jac->printstatistics) {
280390e7148SBarry Smith     args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1";
281390e7148SBarry Smith     args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1";
28216d9e3a6SLisandro Dalcin   }
28316d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
28430e6f737SJed Brown   if (cnt) PetscStackCallHypre(0,HYPRE_EuclidSetParams,(jac->hsolver,cnt,args));
28516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
28616d9e3a6SLisandro Dalcin }
28716d9e3a6SLisandro Dalcin 
28816d9e3a6SLisandro Dalcin #undef __FUNCT__
28916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Euclid"
29016d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
29116d9e3a6SLisandro Dalcin {
29216d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
29316d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
294ace3abfcSBarry Smith   PetscBool      iascii;
29516d9e3a6SLisandro Dalcin 
29616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
297251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
29816d9e3a6SLisandro Dalcin   if (iascii) {
29916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n");CHKERRQ(ierr);
30016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr);
30116d9e3a6SLisandro Dalcin     if (jac->bjilu) {
30216d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr);
30316d9e3a6SLisandro Dalcin     }
30416d9e3a6SLisandro Dalcin   }
30516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
30616d9e3a6SLisandro Dalcin }
30716d9e3a6SLisandro Dalcin 
30816d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
30916d9e3a6SLisandro Dalcin 
31016d9e3a6SLisandro Dalcin #undef __FUNCT__
31116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG"
31216d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
31316d9e3a6SLisandro Dalcin {
31416d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
31516d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
31616d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
31716d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
31816d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
31916d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
3204ddd07fcSJed Brown   PetscInt           hierr;
32116d9e3a6SLisandro Dalcin 
32216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
32316d9e3a6SLisandro Dalcin   ierr = VecSet(x,0.0);CHKERRQ(ierr);
32416d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
32516d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
32616d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
32716d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
32816d9e3a6SLisandro Dalcin 
32930e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
33030e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJVectorGetObject,(jac->b,(void**)&jbv));
33130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJVectorGetObject,(jac->x,(void**)&jxv));
33216d9e3a6SLisandro Dalcin 
33316d9e3a6SLisandro Dalcin   hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv);
33416d9e3a6SLisandro Dalcin   /* error code of 1 in BoomerAMG merely means convergence not achieved */
335e32f2f54SBarry Smith   if (hierr && (hierr != 1)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
33616d9e3a6SLisandro Dalcin   if (hierr) hypre__global_error = 0;
33716d9e3a6SLisandro Dalcin 
33816d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
33916d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
34016d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
34116d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
34216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
34316d9e3a6SLisandro Dalcin }
34416d9e3a6SLisandro Dalcin 
34516d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
3460f1074feSSatish Balay static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"};
34716d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
348*65de4495SJed Brown /* The following corresponds to HYPRE_BoomerAMGSetRelaxType which has many missing numbers in the enum */
349*65de4495SJed Brown static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","seqboundary-Gauss-Seidel","SOR/Jacobi","backward-SOR/Jacobi",
350*65de4495SJed Brown                                                   "" /* [5] hybrid chaotic Gauss-Seidel (works only with OpenMP) */,"symmetric-SOR/Jacobi",
351*65de4495SJed Brown                                                   "" /* 7 */,"l1scaled-SOR/Jacobi","Gaussian-elimination",
352*65de4495SJed Brown                                                   "" /* 10 */, "" /* 11 */, "" /* 12 */, "" /* 13 */, "" /* 14 */,
353*65de4495SJed Brown                                                   "CG" /* non-stationary */,"Chebyshev","FCF-Jacobi","l1scaled-Jacobi"};
3540f1074feSSatish Balay static const char *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i",
3550f1074feSSatish Balay                                                   "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"};
35616d9e3a6SLisandro Dalcin #undef __FUNCT__
35716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG"
35816d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
35916d9e3a6SLisandro Dalcin {
36016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
36116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
3624ddd07fcSJed Brown   PetscInt       n,indx,level;
363ace3abfcSBarry Smith   PetscBool      flg, tmp_truth;
36416d9e3a6SLisandro Dalcin   double         tmpdbl, twodbl[2];
36516d9e3a6SLisandro Dalcin 
36616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
36716d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr);
3684336a9eeSBarry Smith   ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType+1,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr);
36916d9e3a6SLisandro Dalcin   if (flg) {
3704336a9eeSBarry Smith     jac->cycletype = indx+1;
37130e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
37216d9e3a6SLisandro Dalcin   }
37316d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr);
37416d9e3a6SLisandro Dalcin   if (flg) {
37565e19b50SBarry Smith     if (jac->maxlevels < 2) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
37630e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
37716d9e3a6SLisandro Dalcin   }
37816d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr);
37916d9e3a6SLisandro Dalcin   if (flg) {
38065e19b50SBarry Smith     if (jac->maxiter < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
38130e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
38216d9e3a6SLisandro Dalcin   }
3830f1074feSSatish 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);
38416d9e3a6SLisandro Dalcin   if (flg) {
38565e19b50SBarry 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);
38630e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
38716d9e3a6SLisandro Dalcin   }
38816d9e3a6SLisandro Dalcin 
3890f1074feSSatish Balay   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr);
39016d9e3a6SLisandro Dalcin   if (flg) {
39165e19b50SBarry 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);
39230e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
39316d9e3a6SLisandro Dalcin   }
39416d9e3a6SLisandro Dalcin 
3950f1074feSSatish 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);
3960f1074feSSatish Balay   if (flg) {
39765e19b50SBarry 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);
39830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
3990f1074feSSatish Balay   }
4000f1074feSSatish Balay 
4010f1074feSSatish Balay  ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);CHKERRQ(ierr);
4020f1074feSSatish Balay   if (flg) {
40365e19b50SBarry 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);
4040f1074feSSatish Balay 
40530e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
4060f1074feSSatish Balay   }
4070f1074feSSatish Balay 
4080f1074feSSatish Balay 
4090f1074feSSatish 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);
4100f1074feSSatish Balay   if (flg) {
41165e19b50SBarry 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);
4120f1074feSSatish Balay 
41330e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
4140f1074feSSatish Balay   }
4150f1074feSSatish Balay 
4160f1074feSSatish Balay 
41716d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr);
41816d9e3a6SLisandro Dalcin   if (flg) {
41965e19b50SBarry 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);
42030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
42116d9e3a6SLisandro Dalcin   }
42216d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr);
42316d9e3a6SLisandro Dalcin   if (flg) {
42465e19b50SBarry 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);
42565e19b50SBarry 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);
42630e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
42716d9e3a6SLisandro Dalcin   }
42816d9e3a6SLisandro Dalcin 
42916d9e3a6SLisandro Dalcin   /* Grid sweeps */
4300f1074feSSatish 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);
43116d9e3a6SLisandro Dalcin   if (flg) {
43230e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver,indx));
43316d9e3a6SLisandro Dalcin     /* modify the jac structure so we can view the updated options with PC_View */
43416d9e3a6SLisandro Dalcin     jac->gridsweeps[0] = indx;
4350f1074feSSatish Balay     jac->gridsweeps[1] = indx;
4360f1074feSSatish Balay     /*defaults coarse to 1 */
4370f1074feSSatish Balay     jac->gridsweeps[2] = 1;
43816d9e3a6SLisandro Dalcin   }
4390f1074feSSatish Balay 
4400f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx ,&flg);CHKERRQ(ierr);
44116d9e3a6SLisandro Dalcin   if (flg) {
44230e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 1));
4430f1074feSSatish Balay     jac->gridsweeps[0] = indx;
44416d9e3a6SLisandro Dalcin   }
44516d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr);
44616d9e3a6SLisandro Dalcin   if (flg) {
44730e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 2));
4480f1074feSSatish Balay     jac->gridsweeps[1] = indx;
44916d9e3a6SLisandro Dalcin   }
4500f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr);
45116d9e3a6SLisandro Dalcin   if (flg) {
45230e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 3));
4530f1074feSSatish Balay     jac->gridsweeps[2] = indx;
45416d9e3a6SLisandro Dalcin   }
45516d9e3a6SLisandro Dalcin 
45616d9e3a6SLisandro Dalcin   /* Relax type */
4570f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_all","Relax type for the up and down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
45816d9e3a6SLisandro Dalcin   if (flg) {
4590f1074feSSatish Balay     jac->relaxtype[0] = jac->relaxtype[1]  = indx;
46030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, indx));
4610f1074feSSatish Balay     /* by default, coarse type set to 9 */
4620f1074feSSatish Balay     jac->relaxtype[2] = 9;
4630f1074feSSatish Balay 
46416d9e3a6SLisandro Dalcin   }
4650f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
46616d9e3a6SLisandro Dalcin   if (flg) {
46716d9e3a6SLisandro Dalcin     jac->relaxtype[0] = indx;
46830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 1));
46916d9e3a6SLisandro Dalcin   }
4700f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
47116d9e3a6SLisandro Dalcin   if (flg) {
4720f1074feSSatish Balay     jac->relaxtype[1] = indx;
47330e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 2));
47416d9e3a6SLisandro Dalcin   }
47516d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr);
47616d9e3a6SLisandro Dalcin   if (flg) {
4770f1074feSSatish Balay     jac->relaxtype[2] = indx;
47830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 3));
47916d9e3a6SLisandro Dalcin   }
48016d9e3a6SLisandro Dalcin 
48116d9e3a6SLisandro Dalcin   /* Relaxation Weight */
48216d9e3a6SLisandro 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);
48316d9e3a6SLisandro Dalcin   if (flg) {
48430e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetRelaxWt,(jac->hsolver,tmpdbl));
48516d9e3a6SLisandro Dalcin     jac->relaxweight = tmpdbl;
48616d9e3a6SLisandro Dalcin   }
48716d9e3a6SLisandro Dalcin 
48816d9e3a6SLisandro Dalcin   n=2;
48916d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
49016d9e3a6SLisandro 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);
49116d9e3a6SLisandro Dalcin   if (flg) {
49216d9e3a6SLisandro Dalcin     if (n == 2) {
49316d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
49430e6f737SJed Brown       PetscStackCallHypre(0,HYPRE_BoomerAMGSetLevelRelaxWt,(jac->hsolver,twodbl[0],indx));
49516d9e3a6SLisandro Dalcin     } else {
49665e19b50SBarry 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);
49716d9e3a6SLisandro Dalcin     }
49816d9e3a6SLisandro Dalcin   }
49916d9e3a6SLisandro Dalcin 
50016d9e3a6SLisandro Dalcin   /* Outer relaxation Weight */
50116d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_boomeramg_outer_relax_weight_all","Outer relaxation weight for all levels ( -k = determined with k CG steps)","None",jac->outerrelaxweight, &tmpdbl ,&flg);CHKERRQ(ierr);
50216d9e3a6SLisandro Dalcin   if (flg) {
50330e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetOuterWt,( jac->hsolver, tmpdbl));
50416d9e3a6SLisandro Dalcin     jac->outerrelaxweight = tmpdbl;
50516d9e3a6SLisandro Dalcin   }
50616d9e3a6SLisandro Dalcin 
50716d9e3a6SLisandro Dalcin   n=2;
50816d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
50916d9e3a6SLisandro Dalcin   ierr = PetscOptionsRealArray("-pc_hypre_boomeramg_outer_relax_weight_level","Set the outer relaxation weight for a particular level (weight,level)","None",twodbl, &n, &flg);CHKERRQ(ierr);
51016d9e3a6SLisandro Dalcin   if (flg) {
51116d9e3a6SLisandro Dalcin     if (n == 2) {
51216d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
51330e6f737SJed Brown       PetscStackCallHypre(0,HYPRE_BoomerAMGSetLevelOuterWt,( jac->hsolver, twodbl[0], indx));
51416d9e3a6SLisandro Dalcin     } else {
51565e19b50SBarry 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);
51616d9e3a6SLisandro Dalcin     }
51716d9e3a6SLisandro Dalcin   }
51816d9e3a6SLisandro Dalcin 
51916d9e3a6SLisandro Dalcin   /* the Relax Order */
520acfcf0e5SJed Brown   ierr = PetscOptionsBool( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
52116d9e3a6SLisandro Dalcin 
52216d9e3a6SLisandro Dalcin   if (flg) {
52316d9e3a6SLisandro Dalcin     jac->relaxorder = 0;
52430e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
52516d9e3a6SLisandro Dalcin   }
52616d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr);
52716d9e3a6SLisandro Dalcin   if (flg) {
52816d9e3a6SLisandro Dalcin     jac->measuretype = indx;
52930e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
53016d9e3a6SLisandro Dalcin   }
5310f1074feSSatish Balay   /* update list length 3/07 */
5320f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,11,HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr);
53316d9e3a6SLisandro Dalcin   if (flg) {
53416d9e3a6SLisandro Dalcin     jac->coarsentype = indx;
53530e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
53616d9e3a6SLisandro Dalcin   }
5370f1074feSSatish Balay 
5380f1074feSSatish Balay   /* new 3/07 */
5390f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,14,HYPREBoomerAMGInterpType[0],&indx,&flg);CHKERRQ(ierr);
5400f1074feSSatish Balay   if (flg) {
5410f1074feSSatish Balay     jac->interptype = indx;
54230e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
5430f1074feSSatish Balay   }
5440f1074feSSatish Balay 
545b96a4a96SBarry Smith   ierr = PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);CHKERRQ(ierr);
54616d9e3a6SLisandro Dalcin   if (flg) {
547b96a4a96SBarry Smith     level = 3;
54816d9e3a6SLisandro Dalcin     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
549b96a4a96SBarry Smith     jac->printstatistics = PETSC_TRUE;
55030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetPrintLevel,(jac->hsolver,level));
5512ae77aedSBarry Smith   }
5522ae77aedSBarry Smith 
553b96a4a96SBarry Smith   ierr = PetscOptionsName("-pc_hypre_boomeramg_print_debug","Print debug information","None",&flg);CHKERRQ(ierr);
5542ae77aedSBarry Smith   if (flg) {
555b96a4a96SBarry Smith     level = 3;
5562ae77aedSBarry Smith     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
557b96a4a96SBarry Smith     jac->printstatistics = PETSC_TRUE;
55830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetDebugFlag,(jac->hsolver,level));
55916d9e3a6SLisandro Dalcin   }
5608f87f92bSBarry Smith 
561acfcf0e5SJed Brown   ierr = PetscOptionsBool( "-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5628f87f92bSBarry Smith   if (flg && tmp_truth) {
5638f87f92bSBarry Smith     jac->nodal_coarsen = 1;
56430e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetNodal,(jac->hsolver,1));
5658f87f92bSBarry Smith   }
5668f87f92bSBarry Smith 
567acfcf0e5SJed Brown   ierr = PetscOptionsBool( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5688f87f92bSBarry Smith   if (flg && tmp_truth) {
5698f87f92bSBarry Smith     PetscInt tmp_int;
5708f87f92bSBarry Smith     ierr = PetscOptionsInt( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);CHKERRQ(ierr);
5718f87f92bSBarry Smith     if (flg) jac->nodal_relax_levels = tmp_int;
57230e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetSmoothType,(jac->hsolver,6));
57330e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetDomainType,(jac->hsolver,1));
57430e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetOverlap,(jac->hsolver,0));
57530e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetSmoothNumLevels,(jac->hsolver,jac->nodal_relax_levels));
5768f87f92bSBarry Smith   }
5778f87f92bSBarry Smith 
57816d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
57916d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
58016d9e3a6SLisandro Dalcin }
58116d9e3a6SLisandro Dalcin 
58216d9e3a6SLisandro Dalcin #undef __FUNCT__
58316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG"
584ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool  guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
58516d9e3a6SLisandro Dalcin {
58616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
58716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
5884ddd07fcSJed Brown   PetscInt       oits;
58916d9e3a6SLisandro Dalcin 
59016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
59130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,its*jac->maxiter));
59230e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_BoomerAMGSetTol,(jac->hsolver,rtol));
59316d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_TRUE;
59416d9e3a6SLisandro Dalcin   ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr);
59516d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_FALSE;
59630e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_BoomerAMGGetNumIterations,(jac->hsolver,&oits));
5974d0a8057SBarry Smith   *outits = oits;
5984d0a8057SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
5994d0a8057SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
60030e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
60130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
60216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
60316d9e3a6SLisandro Dalcin }
60416d9e3a6SLisandro Dalcin 
60516d9e3a6SLisandro Dalcin 
60616d9e3a6SLisandro Dalcin #undef __FUNCT__
60716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_BoomerAMG"
60816d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
60916d9e3a6SLisandro Dalcin {
61016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
61116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
612ace3abfcSBarry Smith   PetscBool      iascii;
61316d9e3a6SLisandro Dalcin 
61416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
615251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
61616d9e3a6SLisandro Dalcin   if (iascii) {
61716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr);
61816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr);
61916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr);
62016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr);
62116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr);
62216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr);
6230f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation truncation factor %G\n",jac->truncfactor);CHKERRQ(ierr);
6240f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);CHKERRQ(ierr);
6250f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);CHKERRQ(ierr);
6260f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);CHKERRQ(ierr);
6270f1074feSSatish Balay 
62816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr);
62916d9e3a6SLisandro Dalcin 
6300f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps down         %d\n",jac->gridsweeps[0]);CHKERRQ(ierr);
6310f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps up           %d\n",jac->gridsweeps[1]);CHKERRQ(ierr);
6320f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on coarse    %d\n",jac->gridsweeps[2]);CHKERRQ(ierr);
63316d9e3a6SLisandro Dalcin 
6340f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr);
6350f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr);
6360f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr);
63716d9e3a6SLisandro Dalcin 
63816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax weight  (all)      %G\n",jac->relaxweight);CHKERRQ(ierr);
63916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr);
64016d9e3a6SLisandro Dalcin 
64116d9e3a6SLisandro Dalcin     if (jac->relaxorder) {
64216d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr);
64316d9e3a6SLisandro Dalcin     } else {
64416d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr);
64516d9e3a6SLisandro Dalcin     }
64616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr);
64716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr);
6480f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation type  %s\n",HYPREBoomerAMGInterpType[jac->interptype]);CHKERRQ(ierr);
6498f87f92bSBarry Smith     if (jac->nodal_coarsen) {
6508f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal())\n");CHKERRQ(ierr);
6518f87f92bSBarry Smith     }
6528f87f92bSBarry Smith     if (jac->nodal_relax) {
6538f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);CHKERRQ(ierr);
6548f87f92bSBarry Smith     }
65516d9e3a6SLisandro Dalcin   }
65616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
65716d9e3a6SLisandro Dalcin }
65816d9e3a6SLisandro Dalcin 
65916d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
66016d9e3a6SLisandro Dalcin #undef __FUNCT__
66116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails"
66216d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc)
66316d9e3a6SLisandro Dalcin {
66416d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
66516d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
6664ddd07fcSJed Brown   PetscInt       indx;
667ace3abfcSBarry Smith   PetscBool      flag;
66816d9e3a6SLisandro Dalcin   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};
66916d9e3a6SLisandro Dalcin 
67016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
67116d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr);
67216d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr);
67316d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr);
67416d9e3a6SLisandro Dalcin   if (flag) {
67530e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels));
67616d9e3a6SLisandro Dalcin   }
67716d9e3a6SLisandro Dalcin 
67816d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr);
67916d9e3a6SLisandro Dalcin   if (flag) {
68030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter));
68116d9e3a6SLisandro Dalcin   }
68216d9e3a6SLisandro Dalcin 
68316d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr);
68416d9e3a6SLisandro Dalcin   if (flag) {
68530e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal));
68616d9e3a6SLisandro Dalcin   }
68716d9e3a6SLisandro Dalcin 
688acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_parasails_logging","Print info to screen","None",(PetscBool)jac->logging,(PetscBool*)&jac->logging,&flag);CHKERRQ(ierr);
68916d9e3a6SLisandro Dalcin   if (flag) {
69030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging));
69116d9e3a6SLisandro Dalcin   }
69216d9e3a6SLisandro Dalcin 
693acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscBool)jac->ruse,(PetscBool*)&jac->ruse,&flag);CHKERRQ(ierr);
69416d9e3a6SLisandro Dalcin   if (flag) {
69530e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse));
69616d9e3a6SLisandro Dalcin   }
69716d9e3a6SLisandro Dalcin 
69816d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);CHKERRQ(ierr);
69916d9e3a6SLisandro Dalcin   if (flag) {
70016d9e3a6SLisandro Dalcin     jac->symt = indx;
70130e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
70216d9e3a6SLisandro Dalcin   }
70316d9e3a6SLisandro Dalcin 
70416d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
70516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
70616d9e3a6SLisandro Dalcin }
70716d9e3a6SLisandro Dalcin 
70816d9e3a6SLisandro Dalcin #undef __FUNCT__
70916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_ParaSails"
71016d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
71116d9e3a6SLisandro Dalcin {
71216d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
71316d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
714ace3abfcSBarry Smith   PetscBool      iascii;
71516d9e3a6SLisandro Dalcin   const char     *symt = 0;;
71616d9e3a6SLisandro Dalcin 
71716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
718251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
71916d9e3a6SLisandro Dalcin   if (iascii) {
72016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");CHKERRQ(ierr);
72116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr);
72216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr);
72316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr);
72416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr);
725ace3abfcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscBools[jac->ruse]);CHKERRQ(ierr);
726ace3abfcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscBools[jac->logging]);CHKERRQ(ierr);
72716d9e3a6SLisandro Dalcin     if (!jac->symt) {
72816d9e3a6SLisandro Dalcin       symt = "nonsymmetric matrix and preconditioner";
72916d9e3a6SLisandro Dalcin     } else if (jac->symt == 1) {
73016d9e3a6SLisandro Dalcin       symt = "SPD matrix and preconditioner";
73116d9e3a6SLisandro Dalcin     } else if (jac->symt == 2) {
73216d9e3a6SLisandro Dalcin       symt = "nonsymmetric matrix but SPD preconditioner";
73316d9e3a6SLisandro Dalcin     } else {
73465e19b50SBarry Smith       SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
73516d9e3a6SLisandro Dalcin     }
73616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr);
73716d9e3a6SLisandro Dalcin   }
73816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
73916d9e3a6SLisandro Dalcin }
74016d9e3a6SLisandro Dalcin /* ---------------------------------------------------------------------------------*/
74116d9e3a6SLisandro Dalcin 
74216d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
74316d9e3a6SLisandro Dalcin #undef __FUNCT__
74416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType_HYPRE"
7457087cfbeSBarry Smith PetscErrorCode  PCHYPREGetType_HYPRE(PC pc,const char *name[])
74616d9e3a6SLisandro Dalcin {
74716d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
74816d9e3a6SLisandro Dalcin 
74916d9e3a6SLisandro Dalcin   PetscFunctionBegin;
75016d9e3a6SLisandro Dalcin   *name = jac->hypre_type;
75116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
75216d9e3a6SLisandro Dalcin }
75316d9e3a6SLisandro Dalcin EXTERN_C_END
75416d9e3a6SLisandro Dalcin 
75516d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
75616d9e3a6SLisandro Dalcin #undef __FUNCT__
75716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType_HYPRE"
7587087cfbeSBarry Smith PetscErrorCode  PCHYPRESetType_HYPRE(PC pc,const char name[])
75916d9e3a6SLisandro Dalcin {
76016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
76116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
762ace3abfcSBarry Smith   PetscBool      flag;
76316d9e3a6SLisandro Dalcin 
76416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
76516d9e3a6SLisandro Dalcin   if (jac->hypre_type) {
76616d9e3a6SLisandro Dalcin     ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr);
767e7e72b3dSBarry Smith     if (!flag) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
76816d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
76916d9e3a6SLisandro Dalcin   } else {
77016d9e3a6SLisandro Dalcin     ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr);
77116d9e3a6SLisandro Dalcin   }
77216d9e3a6SLisandro Dalcin 
77316d9e3a6SLisandro Dalcin   jac->maxiter            = PETSC_DEFAULT;
77416d9e3a6SLisandro Dalcin   jac->tol                = PETSC_DEFAULT;
77516d9e3a6SLisandro Dalcin   jac->printstatistics    = PetscLogPrintInfo;
77616d9e3a6SLisandro Dalcin 
77716d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr);
77816d9e3a6SLisandro Dalcin   if (flag) {
77930e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParCSRPilutCreate,(jac->comm_hypre,&jac->hsolver));
78016d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
78116d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Pilut;
78216d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParCSRPilutDestroy;
78316d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParCSRPilutSetup;
78416d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParCSRPilutSolve;
78516d9e3a6SLisandro Dalcin     jac->factorrowsize      = PETSC_DEFAULT;
78616d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
78716d9e3a6SLisandro Dalcin   }
78816d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr);
78916d9e3a6SLisandro Dalcin   if (flag) {
79030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsCreate,(jac->comm_hypre,&jac->hsolver));
79116d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
79216d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_ParaSails;
79316d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParaSailsDestroy;
79416d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParaSailsSetup;
79516d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParaSailsSolve;
79616d9e3a6SLisandro Dalcin     /* initialize */
79716d9e3a6SLisandro Dalcin     jac->nlevels     = 1;
79816d9e3a6SLisandro Dalcin     jac->threshhold  = .1;
79916d9e3a6SLisandro Dalcin     jac->filter      = .1;
80016d9e3a6SLisandro Dalcin     jac->loadbal     = 0;
80116d9e3a6SLisandro Dalcin     if (PetscLogPrintInfo) {
80216d9e3a6SLisandro Dalcin       jac->logging     = (int) PETSC_TRUE;
80316d9e3a6SLisandro Dalcin     } else {
80416d9e3a6SLisandro Dalcin       jac->logging     = (int) PETSC_FALSE;
80516d9e3a6SLisandro Dalcin     }
80616d9e3a6SLisandro Dalcin     jac->ruse = (int) PETSC_FALSE;
80716d9e3a6SLisandro Dalcin     jac->symt   = 0;
80830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels));
80930e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter));
81030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal));
81130e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging));
81230e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse));
81330e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
81416d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
81516d9e3a6SLisandro Dalcin   }
81616d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr);
81716d9e3a6SLisandro Dalcin   if (flag) {
81816d9e3a6SLisandro Dalcin     ierr                    = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver);
81916d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
82016d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Euclid;
82116d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_EuclidDestroy;
82216d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_EuclidSetup;
82316d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_EuclidSolve;
82416d9e3a6SLisandro Dalcin     /* initialization */
82516d9e3a6SLisandro Dalcin     jac->bjilu              = PETSC_FALSE;
82616d9e3a6SLisandro Dalcin     jac->levels             = 1;
82716d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
82816d9e3a6SLisandro Dalcin   }
82916d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr);
83016d9e3a6SLisandro Dalcin   if (flag) {
83116d9e3a6SLisandro Dalcin     ierr                     = HYPRE_BoomerAMGCreate(&jac->hsolver);
83216d9e3a6SLisandro Dalcin     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
83316d9e3a6SLisandro Dalcin     pc->ops->view            = PCView_HYPRE_BoomerAMG;
83416d9e3a6SLisandro Dalcin     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
83516d9e3a6SLisandro Dalcin     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
83616d9e3a6SLisandro Dalcin     jac->destroy             = HYPRE_BoomerAMGDestroy;
83716d9e3a6SLisandro Dalcin     jac->setup               = HYPRE_BoomerAMGSetup;
83816d9e3a6SLisandro Dalcin     jac->solve               = HYPRE_BoomerAMGSolve;
83916d9e3a6SLisandro Dalcin     jac->applyrichardson     = PETSC_FALSE;
84016d9e3a6SLisandro Dalcin     /* these defaults match the hypre defaults */
84116d9e3a6SLisandro Dalcin     jac->cycletype        = 1;
84216d9e3a6SLisandro Dalcin     jac->maxlevels        = 25;
84316d9e3a6SLisandro Dalcin     jac->maxiter          = 1;
8448f87f92bSBarry Smith     jac->tol              = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */
84516d9e3a6SLisandro Dalcin     jac->truncfactor      = 0.0;
84616d9e3a6SLisandro Dalcin     jac->strongthreshold  = .25;
84716d9e3a6SLisandro Dalcin     jac->maxrowsum        = .9;
84816d9e3a6SLisandro Dalcin     jac->coarsentype      = 6;
84916d9e3a6SLisandro Dalcin     jac->measuretype      = 0;
8500f1074feSSatish Balay     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = 1;
8518f87f92bSBarry Smith     jac->relaxtype[0]     = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */
8520f1074feSSatish Balay     jac->relaxtype[2]     = 9; /*G.E. */
85316d9e3a6SLisandro Dalcin     jac->relaxweight      = 1.0;
85416d9e3a6SLisandro Dalcin     jac->outerrelaxweight = 1.0;
85516d9e3a6SLisandro Dalcin     jac->relaxorder       = 1;
8560f1074feSSatish Balay     jac->interptype       = 0;
8570f1074feSSatish Balay     jac->agg_nl           = 0;
8580f1074feSSatish Balay     jac->pmax             = 0;
8590f1074feSSatish Balay     jac->truncfactor      = 0.0;
8600f1074feSSatish Balay     jac->agg_num_paths    = 1;
8618f87f92bSBarry Smith 
8628f87f92bSBarry Smith     jac->nodal_coarsen    = 0;
8638f87f92bSBarry Smith     jac->nodal_relax      = PETSC_FALSE;
8648f87f92bSBarry Smith     jac->nodal_relax_levels = 1;
86530e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
86630e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
86730e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
86830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
86930e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
87030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
87130e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
87230e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
87330e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
87430e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
87530e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
87630e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
87730e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
87830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
87930e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, jac->relaxtype[0]));  /*defaults coarse to 9*/
88030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver, jac->gridsweeps[0])); /*defaults coarse to 1 */
88116d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
88216d9e3a6SLisandro Dalcin   }
883503cfb0cSBarry Smith   ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr);
88416d9e3a6SLisandro Dalcin   jac->hypre_type = PETSC_NULL;
88565e19b50SBarry Smith   SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
88616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
88716d9e3a6SLisandro Dalcin }
88816d9e3a6SLisandro Dalcin EXTERN_C_END
88916d9e3a6SLisandro Dalcin 
89016d9e3a6SLisandro Dalcin /*
89116d9e3a6SLisandro Dalcin     It only gets here if the HYPRE type has not been set before the call to
89216d9e3a6SLisandro Dalcin    ...SetFromOptions() which actually is most of the time
89316d9e3a6SLisandro Dalcin */
89416d9e3a6SLisandro Dalcin #undef __FUNCT__
89516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE"
89616d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE(PC pc)
89716d9e3a6SLisandro Dalcin {
89816d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
8994ddd07fcSJed Brown   PetscInt       indx;
90016d9e3a6SLisandro Dalcin   const char     *type[] = {"pilut","parasails","boomeramg","euclid"};
901ace3abfcSBarry Smith   PetscBool      flg;
90216d9e3a6SLisandro Dalcin 
90316d9e3a6SLisandro Dalcin   PetscFunctionBegin;
90416d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr);
90502a17cd4SBarry Smith   ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"boomeramg",&indx,&flg);CHKERRQ(ierr);
90616d9e3a6SLisandro Dalcin   if (flg) {
90716d9e3a6SLisandro Dalcin     ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr);
90802a17cd4SBarry Smith   } else {
90902a17cd4SBarry Smith     ierr = PCHYPRESetType_HYPRE(pc,"boomeramg");CHKERRQ(ierr);
91016d9e3a6SLisandro Dalcin   }
91116d9e3a6SLisandro Dalcin   if (pc->ops->setfromoptions) {
91216d9e3a6SLisandro Dalcin     ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr);
91316d9e3a6SLisandro Dalcin   }
91416d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
91516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
91616d9e3a6SLisandro Dalcin }
91716d9e3a6SLisandro Dalcin 
91816d9e3a6SLisandro Dalcin #undef __FUNCT__
91916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType"
92016d9e3a6SLisandro Dalcin /*@C
92116d9e3a6SLisandro Dalcin      PCHYPRESetType - Sets which hypre preconditioner you wish to use
92216d9e3a6SLisandro Dalcin 
92316d9e3a6SLisandro Dalcin    Input Parameters:
92416d9e3a6SLisandro Dalcin +     pc - the preconditioner context
92516d9e3a6SLisandro Dalcin -     name - either  pilut, parasails, boomeramg, euclid
92616d9e3a6SLisandro Dalcin 
92716d9e3a6SLisandro Dalcin    Options Database Keys:
92816d9e3a6SLisandro Dalcin    -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
92916d9e3a6SLisandro Dalcin 
93016d9e3a6SLisandro Dalcin    Level: intermediate
93116d9e3a6SLisandro Dalcin 
93216d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
93316d9e3a6SLisandro Dalcin            PCHYPRE
93416d9e3a6SLisandro Dalcin 
93516d9e3a6SLisandro Dalcin @*/
9367087cfbeSBarry Smith PetscErrorCode  PCHYPRESetType(PC pc,const char name[])
93716d9e3a6SLisandro Dalcin {
9384ac538c5SBarry Smith   PetscErrorCode ierr;
93916d9e3a6SLisandro Dalcin 
94016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9410700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
94216d9e3a6SLisandro Dalcin   PetscValidCharPointer(name,2);
9434ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCHYPRESetType_C",(PC,const char[]),(pc,name));CHKERRQ(ierr);
94416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
94516d9e3a6SLisandro Dalcin }
94616d9e3a6SLisandro Dalcin 
94716d9e3a6SLisandro Dalcin #undef __FUNCT__
94816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType"
94916d9e3a6SLisandro Dalcin /*@C
95016d9e3a6SLisandro Dalcin      PCHYPREGetType - Gets which hypre preconditioner you are using
95116d9e3a6SLisandro Dalcin 
95216d9e3a6SLisandro Dalcin    Input Parameter:
95316d9e3a6SLisandro Dalcin .     pc - the preconditioner context
95416d9e3a6SLisandro Dalcin 
95516d9e3a6SLisandro Dalcin    Output Parameter:
95616d9e3a6SLisandro Dalcin .     name - either  pilut, parasails, boomeramg, euclid
95716d9e3a6SLisandro Dalcin 
95816d9e3a6SLisandro Dalcin    Level: intermediate
95916d9e3a6SLisandro Dalcin 
96016d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC,
96116d9e3a6SLisandro Dalcin            PCHYPRE
96216d9e3a6SLisandro Dalcin 
96316d9e3a6SLisandro Dalcin @*/
9647087cfbeSBarry Smith PetscErrorCode  PCHYPREGetType(PC pc,const char *name[])
96516d9e3a6SLisandro Dalcin {
9664ac538c5SBarry Smith   PetscErrorCode ierr;
96716d9e3a6SLisandro Dalcin 
96816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9690700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
97016d9e3a6SLisandro Dalcin   PetscValidPointer(name,2);
9714ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCHYPREGetType_C",(PC,const char*[]),(pc,name));CHKERRQ(ierr);
97216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
97316d9e3a6SLisandro Dalcin }
97416d9e3a6SLisandro Dalcin 
97516d9e3a6SLisandro Dalcin /*MC
97616d9e3a6SLisandro Dalcin      PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre
97716d9e3a6SLisandro Dalcin 
97816d9e3a6SLisandro Dalcin    Options Database Keys:
97916d9e3a6SLisandro Dalcin +   -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
98016d9e3a6SLisandro Dalcin -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
98116d9e3a6SLisandro Dalcin           preconditioner
98216d9e3a6SLisandro Dalcin 
98316d9e3a6SLisandro Dalcin    Level: intermediate
98416d9e3a6SLisandro Dalcin 
98516d9e3a6SLisandro Dalcin    Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()),
98616d9e3a6SLisandro Dalcin           the many hypre options can ONLY be set via the options database (e.g. the command line
98716d9e3a6SLisandro Dalcin           or with PetscOptionsSetValue(), there are no functions to set them)
98816d9e3a6SLisandro Dalcin 
98916d9e3a6SLisandro Dalcin           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations
9900f1074feSSatish Balay           (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if
9910f1074feSSatish Balay           -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner
9920f1074feSSatish Balay           (-pc_hypre_boomeramg_rtol should be set to 0.0 - the default - to strictly use a fixed number of
9938f87f92bSBarry Smith           iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations
9940f1074feSSatish Balay           and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10
9950f1074feSSatish Balay           then AT MOST twenty V-cycles of boomeramg will be called.
99616d9e3a6SLisandro Dalcin 
9970f1074feSSatish Balay            Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation
9980f1074feSSatish Balay            (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry.
9990f1074feSSatish Balay            Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi.
100016d9e3a6SLisandro Dalcin           If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly
100116d9e3a6SLisandro Dalcin           and use -ksp_max_it to control the number of V-cycles.
100216d9e3a6SLisandro Dalcin           (see the PETSc FAQ.html at the PETSc website under the Documentation tab).
100316d9e3a6SLisandro Dalcin 
100416d9e3a6SLisandro Dalcin 	  2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option
100516d9e3a6SLisandro Dalcin 	  -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L.
100616d9e3a6SLisandro Dalcin 
10079e5bc791SBarry Smith           See PCPFMG for access to the hypre Struct PFMG solver
10089e5bc791SBarry Smith 
100916d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
10109e5bc791SBarry Smith            PCHYPRESetType(), PCPFMG
101116d9e3a6SLisandro Dalcin 
101216d9e3a6SLisandro Dalcin M*/
101316d9e3a6SLisandro Dalcin 
101416d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
101516d9e3a6SLisandro Dalcin #undef __FUNCT__
101616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCCreate_HYPRE"
10177087cfbeSBarry Smith PetscErrorCode  PCCreate_HYPRE(PC pc)
101816d9e3a6SLisandro Dalcin {
101916d9e3a6SLisandro Dalcin   PC_HYPRE       *jac;
102016d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
102116d9e3a6SLisandro Dalcin 
102216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
102338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_HYPRE,&jac);CHKERRQ(ierr);
102416d9e3a6SLisandro Dalcin   pc->data                 = jac;
102516d9e3a6SLisandro Dalcin   pc->ops->destroy         = PCDestroy_HYPRE;
102616d9e3a6SLisandro Dalcin   pc->ops->setfromoptions  = PCSetFromOptions_HYPRE;
102716d9e3a6SLisandro Dalcin   pc->ops->setup           = PCSetUp_HYPRE;
102816d9e3a6SLisandro Dalcin   pc->ops->apply           = PCApply_HYPRE;
102916d9e3a6SLisandro Dalcin   jac->comm_hypre          = MPI_COMM_NULL;
103016d9e3a6SLisandro Dalcin   jac->hypre_type          = PETSC_NULL;
103116d9e3a6SLisandro Dalcin   /* duplicate communicator for hypre */
10327adad957SLisandro Dalcin   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(jac->comm_hypre));CHKERRQ(ierr);
103316d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr);
103416d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr);
103516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
103616d9e3a6SLisandro Dalcin }
103716d9e3a6SLisandro Dalcin EXTERN_C_END
1038ebc551c0SBarry Smith 
1039f91d8e95SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------*/
1040f91d8e95SBarry Smith 
1041b862ddfaSBarry Smith /* this include is needed ONLY to allow access to the private data inside the Mat object specific to hypre */
1042b45d2f2cSJed Brown #include <petsc-private/matimpl.h>
1043ebc551c0SBarry Smith 
1044ebc551c0SBarry Smith typedef struct {
104568326731SBarry Smith   MPI_Comm            hcomm;       /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */
1046f91d8e95SBarry Smith   HYPRE_StructSolver  hsolver;
10479e5bc791SBarry Smith 
10489e5bc791SBarry Smith   /* keep copy of PFMG options used so may view them */
10494ddd07fcSJed Brown   PetscInt            its;
10509e5bc791SBarry Smith   double              tol;
10514ddd07fcSJed Brown   PetscInt            relax_type;
10524ddd07fcSJed Brown   PetscInt            rap_type;
10534ddd07fcSJed Brown   PetscInt            num_pre_relax,num_post_relax;
10544ddd07fcSJed Brown   PetscInt            max_levels;
1055ebc551c0SBarry Smith } PC_PFMG;
1056ebc551c0SBarry Smith 
1057ebc551c0SBarry Smith #undef __FUNCT__
1058ebc551c0SBarry Smith #define __FUNCT__ "PCDestroy_PFMG"
1059ebc551c0SBarry Smith PetscErrorCode PCDestroy_PFMG(PC pc)
1060ebc551c0SBarry Smith {
1061ebc551c0SBarry Smith   PetscErrorCode ierr;
1062f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1063ebc551c0SBarry Smith 
1064ebc551c0SBarry Smith   PetscFunctionBegin;
106530e6f737SJed Brown   if (ex->hsolver) {PetscStackCallHypre(0,HYPRE_StructPFMGDestroy,(ex->hsolver));}
1066f91d8e95SBarry Smith   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1067c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1068ebc551c0SBarry Smith   PetscFunctionReturn(0);
1069ebc551c0SBarry Smith }
1070ebc551c0SBarry Smith 
10719e5bc791SBarry Smith static const char *PFMGRelaxType[]   = {"Jacobi","Weighted-Jacobi","symmetric-Red/Black-Gauss-Seidel","Red/Black-Gauss-Seidel"};
10729e5bc791SBarry Smith static const char *PFMGRAPType[]   = {"Galerkin","non-Galerkin"};
10739e5bc791SBarry Smith 
1074ebc551c0SBarry Smith #undef __FUNCT__
1075ebc551c0SBarry Smith #define __FUNCT__ "PCView_PFMG"
1076ebc551c0SBarry Smith PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer)
1077ebc551c0SBarry Smith {
1078ebc551c0SBarry Smith   PetscErrorCode ierr;
1079ace3abfcSBarry Smith   PetscBool      iascii;
1080f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1081ebc551c0SBarry Smith 
1082ebc551c0SBarry Smith   PetscFunctionBegin;
1083251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10849e5bc791SBarry Smith   if (iascii) {
10859e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG preconditioning\n");CHKERRQ(ierr);
10869e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
10879e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
10889e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
10899e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: RAP type %s\n",PFMGRAPType[ex->rap_type]);CHKERRQ(ierr);
10909e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
10913b46a515SGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max levels %d\n",ex->max_levels);CHKERRQ(ierr);
10929e5bc791SBarry Smith   }
1093ebc551c0SBarry Smith   PetscFunctionReturn(0);
1094ebc551c0SBarry Smith }
1095ebc551c0SBarry Smith 
10969e5bc791SBarry Smith 
1097ebc551c0SBarry Smith #undef __FUNCT__
1098ebc551c0SBarry Smith #define __FUNCT__ "PCSetFromOptions_PFMG"
1099ebc551c0SBarry Smith PetscErrorCode PCSetFromOptions_PFMG(PC pc)
1100ebc551c0SBarry Smith {
1101ebc551c0SBarry Smith   PetscErrorCode ierr;
1102f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1103ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
1104ebc551c0SBarry Smith 
1105ebc551c0SBarry Smith   PetscFunctionBegin;
1106ebc551c0SBarry Smith   ierr = PetscOptionsHead("PFMG options");CHKERRQ(ierr);
1107acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
110868326731SBarry Smith   if (flg) {
1109a0324ebeSBarry Smith     int level=3;
111030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_StructPFMGSetPrintLevel,(ex->hsolver,level));
111168326731SBarry Smith   }
11129e5bc791SBarry 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);
111330e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetMaxIter,(ex->hsolver,ex->its));
11149e5bc791SBarry 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);
111530e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetNumPreRelax,(ex->hsolver,ex->num_pre_relax));
11169e5bc791SBarry 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);
111730e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetNumPostRelax,(ex->hsolver,ex->num_post_relax));
11189e5bc791SBarry Smith 
11193b46a515SGlenn Hammond   ierr = PetscOptionsInt("-pc_pfmg_max_levels","Max Levels for MG hierarchy","HYPRE_StructPFMGSetMaxLevels",ex->max_levels,&ex->max_levels,PETSC_NULL);CHKERRQ(ierr);
112030e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetMaxLevels,(ex->hsolver,ex->max_levels));
11213b46a515SGlenn Hammond 
11229e5bc791SBarry Smith   ierr = PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr);
112330e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetTol,(ex->hsolver,ex->tol));
11249e5bc791SBarry Smith   ierr = PetscOptionsEList("-pc_pfmg_relax_type","Relax type for the up and down cycles","HYPRE_StructPFMGSetRelaxType",PFMGRelaxType,4,PFMGRelaxType[ex->relax_type],&ex->relax_type,PETSC_NULL);CHKERRQ(ierr);
112530e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetRelaxType,(ex->hsolver, ex->relax_type));
11269e5bc791SBarry Smith   ierr = PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,2,PFMGRAPType[ex->rap_type],&ex->rap_type,PETSC_NULL);CHKERRQ(ierr);
112730e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetRAPType,(ex->hsolver, ex->rap_type));
1128ebc551c0SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1129ebc551c0SBarry Smith   PetscFunctionReturn(0);
1130ebc551c0SBarry Smith }
1131ebc551c0SBarry Smith 
1132f91d8e95SBarry Smith #undef __FUNCT__
1133f91d8e95SBarry Smith #define __FUNCT__ "PCApply_PFMG"
1134f91d8e95SBarry Smith PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y)
1135f91d8e95SBarry Smith {
1136f91d8e95SBarry Smith   PetscErrorCode  ierr;
1137f91d8e95SBarry Smith   PC_PFMG         *ex = (PC_PFMG*) pc->data;
1138f91d8e95SBarry Smith   PetscScalar     *xx,*yy;
11394ddd07fcSJed Brown   PetscInt        ilower[3],iupper[3];
114068326731SBarry Smith   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
1141f91d8e95SBarry Smith 
1142f91d8e95SBarry Smith   PetscFunctionBegin;
1143aa219208SBarry Smith   ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1144f91d8e95SBarry Smith   iupper[0] += ilower[0] - 1;
1145f91d8e95SBarry Smith   iupper[1] += ilower[1] - 1;
1146f91d8e95SBarry Smith   iupper[2] += ilower[2] - 1;
1147f91d8e95SBarry Smith 
1148f91d8e95SBarry Smith   /* copy x values over to hypre */
114930e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructVectorSetConstantValues,(mx->hb,0.0));
1150f91d8e95SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
115130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructVectorSetBoxValues,(mx->hb,ilower,iupper,xx));
1152f91d8e95SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
115330e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructVectorAssemble,(mx->hb));
115430e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSolve,(ex->hsolver,mx->hmat,mx->hb,mx->hx));
1155f91d8e95SBarry Smith 
1156f91d8e95SBarry Smith   /* copy solution values back to PETSc */
1157f91d8e95SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
115830e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructVectorGetBoxValues,(mx->hx,ilower,iupper,yy));
1159f91d8e95SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1160f91d8e95SBarry Smith   PetscFunctionReturn(0);
1161f91d8e95SBarry Smith }
1162f91d8e95SBarry Smith 
11639e5bc791SBarry Smith #undef __FUNCT__
11649e5bc791SBarry Smith #define __FUNCT__ "PCApplyRichardson_PFMG"
1165ace3abfcSBarry 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)
11669e5bc791SBarry Smith {
11679e5bc791SBarry Smith   PC_PFMG        *jac = (PC_PFMG*)pc->data;
11689e5bc791SBarry Smith   PetscErrorCode ierr;
11694ddd07fcSJed Brown   PetscInt       oits;
11709e5bc791SBarry Smith 
11719e5bc791SBarry Smith   PetscFunctionBegin;
117230e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetMaxIter,(jac->hsolver,its*jac->its));
117330e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetTol,(jac->hsolver,rtol));
11749e5bc791SBarry Smith 
11759e5bc791SBarry Smith   ierr = PCApply_PFMG(pc,b,y);CHKERRQ(ierr);
117630e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGGetNumIterations,(jac->hsolver,&oits));
11779e5bc791SBarry Smith   *outits = oits;
11789e5bc791SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
11799e5bc791SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
118030e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetTol,(jac->hsolver,jac->tol));
118130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetMaxIter,(jac->hsolver,jac->its));
11829e5bc791SBarry Smith   PetscFunctionReturn(0);
11839e5bc791SBarry Smith }
11849e5bc791SBarry Smith 
11859e5bc791SBarry Smith 
11863a32d3dbSGlenn Hammond #undef __FUNCT__
11873a32d3dbSGlenn Hammond #define __FUNCT__ "PCSetUp_PFMG"
11883a32d3dbSGlenn Hammond PetscErrorCode PCSetUp_PFMG(PC pc)
11893a32d3dbSGlenn Hammond {
11903a32d3dbSGlenn Hammond   PetscErrorCode  ierr;
11913a32d3dbSGlenn Hammond   PC_PFMG         *ex = (PC_PFMG*) pc->data;
11923a32d3dbSGlenn Hammond   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
1193ace3abfcSBarry Smith   PetscBool       flg;
11943a32d3dbSGlenn Hammond 
11953a32d3dbSGlenn Hammond   PetscFunctionBegin;
1196251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESTRUCT,&flg);CHKERRQ(ierr);
1197e7e72b3dSBarry Smith   if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESTRUCT with this preconditioner");
11983a32d3dbSGlenn Hammond 
11993a32d3dbSGlenn Hammond   /* create the hypre solver object and set its information */
12003a32d3dbSGlenn Hammond   if (ex->hsolver) {
120130e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_StructPFMGDestroy,(ex->hsolver));
12023a32d3dbSGlenn Hammond   }
120330e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
12043a32d3dbSGlenn Hammond   ierr = PCSetFromOptions_PFMG(pc);CHKERRQ(ierr);
120530e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetup,(ex->hsolver,mx->hmat,mx->hb,mx->hx));
120630e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetZeroGuess,(ex->hsolver));
12073a32d3dbSGlenn Hammond   PetscFunctionReturn(0);
12083a32d3dbSGlenn Hammond }
12093a32d3dbSGlenn Hammond 
1210ebc551c0SBarry Smith 
1211ebc551c0SBarry Smith /*MC
1212ebc551c0SBarry Smith      PCPFMG - the hypre PFMG multigrid solver
1213ebc551c0SBarry Smith 
1214ebc551c0SBarry Smith    Level: advanced
1215ebc551c0SBarry Smith 
12169e5bc791SBarry Smith    Options Database:
12179e5bc791SBarry Smith + -pc_pfmg_its <its> number of iterations of PFMG to use as preconditioner
12189e5bc791SBarry Smith . -pc_pfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
12199e5bc791SBarry Smith . -pc_pfmg_num_post_relax <steps> number of smoothing steps after coarse grid
12209e5bc791SBarry Smith . -pc_pfmg_tol <tol> tolerance of PFMG
12219e5bc791SBarry 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
12229e5bc791SBarry Smith - -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin
1223f91d8e95SBarry Smith 
12249e5bc791SBarry Smith    Notes:  This is for CELL-centered descretizations
12259e5bc791SBarry Smith 
12268e395302SJed Brown            This must be used with the MATHYPRESTRUCT matrix type.
1227aa219208SBarry Smith            This is less general than in hypre, it supports only one block per process defined by a PETSc DMDA.
12289e5bc791SBarry Smith 
12299e5bc791SBarry Smith .seealso:  PCMG, MATHYPRESTRUCT
1230ebc551c0SBarry Smith M*/
1231ebc551c0SBarry Smith 
1232ebc551c0SBarry Smith EXTERN_C_BEGIN
1233ebc551c0SBarry Smith #undef __FUNCT__
1234ebc551c0SBarry Smith #define __FUNCT__ "PCCreate_PFMG"
12357087cfbeSBarry Smith PetscErrorCode  PCCreate_PFMG(PC pc)
1236ebc551c0SBarry Smith {
1237ebc551c0SBarry Smith   PetscErrorCode ierr;
1238ebc551c0SBarry Smith   PC_PFMG        *ex;
1239ebc551c0SBarry Smith 
1240ebc551c0SBarry Smith   PetscFunctionBegin;
1241ebc551c0SBarry Smith   ierr = PetscNew(PC_PFMG,&ex);CHKERRQ(ierr);\
124268326731SBarry Smith   pc->data = ex;
1243ebc551c0SBarry Smith 
12449e5bc791SBarry Smith   ex->its            = 1;
12459e5bc791SBarry Smith   ex->tol            = 1.e-8;
12469e5bc791SBarry Smith   ex->relax_type     = 1;
12479e5bc791SBarry Smith   ex->rap_type       = 0;
12489e5bc791SBarry Smith   ex->num_pre_relax  = 1;
12499e5bc791SBarry Smith   ex->num_post_relax = 1;
12503b46a515SGlenn Hammond   ex->max_levels     = 0;
12519e5bc791SBarry Smith 
1252ebc551c0SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_PFMG;
1253ebc551c0SBarry Smith   pc->ops->view            = PCView_PFMG;
1254ebc551c0SBarry Smith   pc->ops->destroy         = PCDestroy_PFMG;
1255f91d8e95SBarry Smith   pc->ops->apply           = PCApply_PFMG;
12569e5bc791SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_PFMG;
125768326731SBarry Smith   pc->ops->setup           = PCSetUp_PFMG;
1258f91d8e95SBarry Smith   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr);
125930e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
1260ebc551c0SBarry Smith   PetscFunctionReturn(0);
1261ebc551c0SBarry Smith }
1262ebc551c0SBarry Smith EXTERN_C_END
1263d851a50bSGlenn Hammond 
1264325fc9f4SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------------------------*/
1265325fc9f4SBarry Smith 
1266d851a50bSGlenn Hammond /* we know we are working with a HYPRE_SStructMatrix */
1267d851a50bSGlenn Hammond typedef struct {
1268d851a50bSGlenn Hammond   MPI_Comm            hcomm;       /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */
1269d851a50bSGlenn Hammond   HYPRE_SStructSolver  ss_solver;
1270d851a50bSGlenn Hammond 
1271d851a50bSGlenn Hammond   /* keep copy of SYSPFMG options used so may view them */
12724ddd07fcSJed Brown   PetscInt            its;
1273d851a50bSGlenn Hammond   double              tol;
12744ddd07fcSJed Brown   PetscInt            relax_type;
12754ddd07fcSJed Brown   PetscInt            num_pre_relax,num_post_relax;
1276d851a50bSGlenn Hammond } PC_SysPFMG;
1277d851a50bSGlenn Hammond 
1278d851a50bSGlenn Hammond #undef __FUNCT__
1279d851a50bSGlenn Hammond #define __FUNCT__ "PCDestroy_SysPFMG"
1280d851a50bSGlenn Hammond PetscErrorCode PCDestroy_SysPFMG(PC pc)
1281d851a50bSGlenn Hammond {
1282d851a50bSGlenn Hammond   PetscErrorCode ierr;
1283d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1284d851a50bSGlenn Hammond 
1285d851a50bSGlenn Hammond   PetscFunctionBegin;
128630e6f737SJed Brown   if (ex->ss_solver) {PetscStackCallHypre(0,HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));}
1287d851a50bSGlenn Hammond   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1288c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1289d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1290d851a50bSGlenn Hammond }
1291d851a50bSGlenn Hammond 
1292d851a50bSGlenn Hammond static const char *SysPFMGRelaxType[]   = {"Weighted-Jacobi","Red/Black-Gauss-Seidel"};
1293d851a50bSGlenn Hammond 
1294d851a50bSGlenn Hammond #undef __FUNCT__
1295d851a50bSGlenn Hammond #define __FUNCT__ "PCView_SysPFMG"
1296d851a50bSGlenn Hammond PetscErrorCode PCView_SysPFMG(PC pc,PetscViewer viewer)
1297d851a50bSGlenn Hammond {
1298d851a50bSGlenn Hammond   PetscErrorCode ierr;
1299ace3abfcSBarry Smith   PetscBool      iascii;
1300d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1301d851a50bSGlenn Hammond 
1302d851a50bSGlenn Hammond   PetscFunctionBegin;
1303251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1304d851a50bSGlenn Hammond   if (iascii) {
1305d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG preconditioning\n");CHKERRQ(ierr);
1306d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
1307d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
1308d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
1309d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
1310d851a50bSGlenn Hammond   }
1311d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1312d851a50bSGlenn Hammond }
1313d851a50bSGlenn Hammond 
1314d851a50bSGlenn Hammond 
1315d851a50bSGlenn Hammond #undef __FUNCT__
1316d851a50bSGlenn Hammond #define __FUNCT__ "PCSetFromOptions_SysPFMG"
1317d851a50bSGlenn Hammond PetscErrorCode PCSetFromOptions_SysPFMG(PC pc)
1318d851a50bSGlenn Hammond {
1319d851a50bSGlenn Hammond   PetscErrorCode ierr;
1320d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1321ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
1322d851a50bSGlenn Hammond 
1323d851a50bSGlenn Hammond   PetscFunctionBegin;
1324d851a50bSGlenn Hammond   ierr = PetscOptionsHead("SysPFMG options");CHKERRQ(ierr);
1325acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_syspfmg_print_statistics","Print statistics","HYPRE_SStructSysPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1326d851a50bSGlenn Hammond   if (flg) {
1327d851a50bSGlenn Hammond     int level=3;
132830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetPrintLevel,(ex->ss_solver,level));
1329d851a50bSGlenn Hammond   }
1330d851a50bSGlenn 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);
133130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetMaxIter,(ex->ss_solver,ex->its));
1332d851a50bSGlenn 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);
133330e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetNumPreRelax,(ex->ss_solver,ex->num_pre_relax));
1334d851a50bSGlenn 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);
133530e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetNumPostRelax,(ex->ss_solver,ex->num_post_relax));
1336d851a50bSGlenn Hammond 
1337d851a50bSGlenn Hammond   ierr = PetscOptionsReal("-pc_syspfmg_tol","Tolerance of SysPFMG","HYPRE_SStructSysPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr);
133830e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetTol,(ex->ss_solver,ex->tol));
1339d851a50bSGlenn 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);
134030e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetRelaxType,(ex->ss_solver, ex->relax_type));
1341d851a50bSGlenn Hammond   ierr = PetscOptionsTail();CHKERRQ(ierr);
1342d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1343d851a50bSGlenn Hammond }
1344d851a50bSGlenn Hammond 
1345d851a50bSGlenn Hammond #undef __FUNCT__
1346d851a50bSGlenn Hammond #define __FUNCT__ "PCApply_SysPFMG"
1347d851a50bSGlenn Hammond PetscErrorCode PCApply_SysPFMG(PC pc,Vec x,Vec y)
1348d851a50bSGlenn Hammond {
1349d851a50bSGlenn Hammond   PetscErrorCode    ierr;
1350d851a50bSGlenn Hammond   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
1351d851a50bSGlenn Hammond   PetscScalar      *xx,*yy;
13524ddd07fcSJed Brown   PetscInt          ilower[3],iupper[3];
1353d851a50bSGlenn Hammond   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct *)(pc->pmat->data);
13544ddd07fcSJed Brown   PetscInt          ordering= mx->dofs_order;
13554ddd07fcSJed Brown   PetscInt          nvars= mx->nvars;
13564ddd07fcSJed Brown   PetscInt          part= 0;
13574ddd07fcSJed Brown   PetscInt          size;
13584ddd07fcSJed Brown   PetscInt          i;
1359d851a50bSGlenn Hammond 
1360d851a50bSGlenn Hammond   PetscFunctionBegin;
1361aa219208SBarry Smith   ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1362d851a50bSGlenn Hammond   iupper[0] += ilower[0] - 1;
1363d851a50bSGlenn Hammond   iupper[1] += ilower[1] - 1;
1364d851a50bSGlenn Hammond   iupper[2] += ilower[2] - 1;
1365d851a50bSGlenn Hammond 
1366d851a50bSGlenn Hammond   size= 1;
1367d851a50bSGlenn Hammond   for (i= 0; i< 3; i++) {
1368d851a50bSGlenn Hammond      size*= (iupper[i]-ilower[i]+1);
1369d851a50bSGlenn Hammond   }
1370d851a50bSGlenn Hammond   /* copy x values over to hypre for variable ordering */
1371d851a50bSGlenn Hammond   if (ordering) {
137230e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
1373d851a50bSGlenn Hammond      ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1374d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
137530e6f737SJed Brown        PetscStackCallHypre(0,HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,ilower,iupper,i,xx+(size*i)));
1376d851a50bSGlenn Hammond      }
1377d851a50bSGlenn Hammond      ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
137830e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_SStructVectorAssemble,(mx->ss_b));
137930e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
138030e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1381d851a50bSGlenn Hammond 
1382d851a50bSGlenn Hammond      /* copy solution values back to PETSc */
1383d851a50bSGlenn Hammond      ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1384d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
138530e6f737SJed Brown        PetscStackCallHypre(0,HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,ilower,iupper,i,yy+(size*i)));
1386d851a50bSGlenn Hammond      }
1387d851a50bSGlenn Hammond      ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1388a65764d7SBarry Smith   } else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
1389d851a50bSGlenn Hammond      PetscScalar     *z;
13904ddd07fcSJed Brown      PetscInt         j, k;
1391d851a50bSGlenn Hammond 
1392d851a50bSGlenn Hammond      ierr = PetscMalloc(nvars*size*sizeof(PetscScalar),&z);CHKERRQ(ierr);
139330e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
1394d851a50bSGlenn Hammond      ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1395d851a50bSGlenn Hammond 
1396d851a50bSGlenn Hammond      /* transform nodal to hypre's variable ordering for sys_pfmg */
1397d851a50bSGlenn Hammond      for (i= 0; i< size; i++) {
1398d851a50bSGlenn Hammond         k= i*nvars;
1399d851a50bSGlenn Hammond         for (j= 0; j< nvars; j++) {
1400d851a50bSGlenn Hammond            z[j*size+i]= xx[k+j];
1401d851a50bSGlenn Hammond         }
1402d851a50bSGlenn Hammond      }
1403d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
140430e6f737SJed Brown        PetscStackCallHypre(0,HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,ilower,iupper,i,z+(size*i)));
1405d851a50bSGlenn Hammond      }
1406d851a50bSGlenn Hammond      ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
140730e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_SStructVectorAssemble,(mx->ss_b));
140830e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1409d851a50bSGlenn Hammond 
1410d851a50bSGlenn Hammond      /* copy solution values back to PETSc */
1411d851a50bSGlenn Hammond      ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1412d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
141330e6f737SJed Brown        PetscStackCallHypre(0,HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,ilower,iupper,i,z+(size*i)));
1414d851a50bSGlenn Hammond      }
1415d851a50bSGlenn Hammond      /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
1416d851a50bSGlenn Hammond      for (i= 0; i< size; i++) {
1417d851a50bSGlenn Hammond         k= i*nvars;
1418d851a50bSGlenn Hammond         for (j= 0; j< nvars; j++) {
1419d851a50bSGlenn Hammond            yy[k+j]= z[j*size+i];
1420d851a50bSGlenn Hammond         }
1421d851a50bSGlenn Hammond      }
1422d851a50bSGlenn Hammond      ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1423d851a50bSGlenn Hammond      ierr = PetscFree(z);CHKERRQ(ierr);
1424d851a50bSGlenn Hammond   }
1425d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1426d851a50bSGlenn Hammond }
1427d851a50bSGlenn Hammond 
1428d851a50bSGlenn Hammond #undef __FUNCT__
1429d851a50bSGlenn Hammond #define __FUNCT__ "PCApplyRichardson_SysPFMG"
1430ace3abfcSBarry 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)
1431d851a50bSGlenn Hammond {
1432d851a50bSGlenn Hammond   PC_SysPFMG    *jac = (PC_SysPFMG*)pc->data;
1433d851a50bSGlenn Hammond   PetscErrorCode ierr;
14344ddd07fcSJed Brown   PetscInt       oits;
1435d851a50bSGlenn Hammond 
1436d851a50bSGlenn Hammond   PetscFunctionBegin;
143730e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,its*jac->its));
143830e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,rtol));
1439d851a50bSGlenn Hammond   ierr = PCApply_SysPFMG(pc,b,y);CHKERRQ(ierr);
144030e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGGetNumIterations,(jac->ss_solver,&oits));
1441d851a50bSGlenn Hammond   *outits = oits;
1442d851a50bSGlenn Hammond   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
1443d851a50bSGlenn Hammond   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
144430e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,jac->tol));
144530e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,jac->its));
1446d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1447d851a50bSGlenn Hammond }
1448d851a50bSGlenn Hammond 
1449d851a50bSGlenn Hammond 
1450d851a50bSGlenn Hammond #undef __FUNCT__
1451d851a50bSGlenn Hammond #define __FUNCT__ "PCSetUp_SysPFMG"
1452d851a50bSGlenn Hammond PetscErrorCode PCSetUp_SysPFMG(PC pc)
1453d851a50bSGlenn Hammond {
1454d851a50bSGlenn Hammond   PetscErrorCode    ierr;
1455d851a50bSGlenn Hammond   PC_SysPFMG        *ex = (PC_SysPFMG*) pc->data;
1456d851a50bSGlenn Hammond   Mat_HYPRESStruct  *mx = (Mat_HYPRESStruct *)(pc->pmat->data);
1457ace3abfcSBarry Smith   PetscBool         flg;
1458d851a50bSGlenn Hammond 
1459d851a50bSGlenn Hammond   PetscFunctionBegin;
1460251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESSTRUCT,&flg);CHKERRQ(ierr);
1461e7e72b3dSBarry Smith   if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESSTRUCT with this preconditioner");
1462d851a50bSGlenn Hammond 
1463d851a50bSGlenn Hammond   /* create the hypre sstruct solver object and set its information */
1464d851a50bSGlenn Hammond   if (ex->ss_solver) {
146530e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));
1466d851a50bSGlenn Hammond   }
146730e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
1468d851a50bSGlenn Hammond   ierr = PCSetFromOptions_SysPFMG(pc);CHKERRQ(ierr);
146930e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetZeroGuess,(ex->ss_solver));
147030e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetup,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1471d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1472d851a50bSGlenn Hammond }
1473d851a50bSGlenn Hammond 
1474d851a50bSGlenn Hammond 
1475d851a50bSGlenn Hammond /*MC
1476d851a50bSGlenn Hammond      PCSysPFMG - the hypre SysPFMG multigrid solver
1477d851a50bSGlenn Hammond 
1478d851a50bSGlenn Hammond    Level: advanced
1479d851a50bSGlenn Hammond 
1480d851a50bSGlenn Hammond    Options Database:
1481d851a50bSGlenn Hammond + -pc_syspfmg_its <its> number of iterations of SysPFMG to use as preconditioner
1482d851a50bSGlenn Hammond . -pc_syspfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
1483d851a50bSGlenn Hammond . -pc_syspfmg_num_post_relax <steps> number of smoothing steps after coarse grid
1484d851a50bSGlenn Hammond . -pc_syspfmg_tol <tol> tolerance of SysPFMG
1485d851a50bSGlenn Hammond . -pc_syspfmg_relax_type -relaxation type for the up and down cycles, one of Weighted-Jacobi,Red/Black-Gauss-Seidel
1486d851a50bSGlenn Hammond 
1487d851a50bSGlenn Hammond    Notes:  This is for CELL-centered descretizations
1488d851a50bSGlenn Hammond 
1489f6680f47SSatish Balay            This must be used with the MATHYPRESSTRUCT matrix type.
1490aa219208SBarry Smith            This is less general than in hypre, it supports only one part, and one block per process defined by a PETSc DMDA.
1491d851a50bSGlenn Hammond            Also, only cell-centered variables.
1492d851a50bSGlenn Hammond 
1493d851a50bSGlenn Hammond .seealso:  PCMG, MATHYPRESSTRUCT
1494d851a50bSGlenn Hammond M*/
1495d851a50bSGlenn Hammond 
1496d851a50bSGlenn Hammond EXTERN_C_BEGIN
1497d851a50bSGlenn Hammond #undef __FUNCT__
1498d851a50bSGlenn Hammond #define __FUNCT__ "PCCreate_SysPFMG"
14997087cfbeSBarry Smith PetscErrorCode  PCCreate_SysPFMG(PC pc)
1500d851a50bSGlenn Hammond {
1501d851a50bSGlenn Hammond   PetscErrorCode     ierr;
1502d851a50bSGlenn Hammond   PC_SysPFMG        *ex;
1503d851a50bSGlenn Hammond 
1504d851a50bSGlenn Hammond   PetscFunctionBegin;
1505d851a50bSGlenn Hammond   ierr = PetscNew(PC_SysPFMG,&ex);CHKERRQ(ierr);\
1506d851a50bSGlenn Hammond   pc->data = ex;
1507d851a50bSGlenn Hammond 
1508d851a50bSGlenn Hammond   ex->its            = 1;
1509d851a50bSGlenn Hammond   ex->tol            = 1.e-8;
1510d851a50bSGlenn Hammond   ex->relax_type     = 1;
1511d851a50bSGlenn Hammond   ex->num_pre_relax  = 1;
1512d851a50bSGlenn Hammond   ex->num_post_relax = 1;
1513d851a50bSGlenn Hammond 
1514d851a50bSGlenn Hammond   pc->ops->setfromoptions  = PCSetFromOptions_SysPFMG;
1515d851a50bSGlenn Hammond   pc->ops->view            = PCView_SysPFMG;
1516d851a50bSGlenn Hammond   pc->ops->destroy         = PCDestroy_SysPFMG;
1517d851a50bSGlenn Hammond   pc->ops->apply           = PCApply_SysPFMG;
1518d851a50bSGlenn Hammond   pc->ops->applyrichardson = PCApplyRichardson_SysPFMG;
1519d851a50bSGlenn Hammond   pc->ops->setup           = PCSetUp_SysPFMG;
1520d851a50bSGlenn Hammond   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr);
152130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
1522d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1523d851a50bSGlenn Hammond }
1524d851a50bSGlenn Hammond EXTERN_C_END
1525