xref: /petsc/src/ksp/pc/impls/hypre/hypre.c (revision 2692d6ee476c3ad7c6573405ba3f5a5bb0473457)
116d9e3a6SLisandro Dalcin #define PETSCKSP_DLL
216d9e3a6SLisandro Dalcin 
316d9e3a6SLisandro Dalcin /*
416d9e3a6SLisandro Dalcin    Provides an interface to the LLNL package hypre
516d9e3a6SLisandro Dalcin */
60f1074feSSatish Balay 
70f1074feSSatish Balay /* Must use hypre 2.0.0 or more recent. */
80f1074feSSatish Balay 
916d9e3a6SLisandro Dalcin #include "private/pcimpl.h"          /*I "petscpc.h" I*/
1016d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
1116d9e3a6SLisandro Dalcin #include "HYPRE.h"
120f1074feSSatish Balay #include "HYPRE_parcsr_ls.h"
130f1074feSSatish Balay #include "_hypre_parcsr_mv.h"
140f1074feSSatish Balay #include "_hypre_IJ_mv.h"
1516d9e3a6SLisandro Dalcin EXTERN_C_END
1616d9e3a6SLisandro Dalcin 
1716d9e3a6SLisandro Dalcin EXTERN PetscErrorCode MatHYPRE_IJMatrixCreate(Mat,HYPRE_IJMatrix*);
1816d9e3a6SLisandro Dalcin EXTERN PetscErrorCode MatHYPRE_IJMatrixCopy(Mat,HYPRE_IJMatrix);
196745373dSBarry Smith EXTERN PetscErrorCode MatHYPRE_IJMatrixFastCopy(Mat,HYPRE_IJMatrix);
2016d9e3a6SLisandro Dalcin EXTERN PetscErrorCode VecHYPRE_IJVectorCreate(Vec,HYPRE_IJVector*);
2116d9e3a6SLisandro Dalcin 
2216d9e3a6SLisandro Dalcin /*
2316d9e3a6SLisandro Dalcin    Private context (data structure) for the  preconditioner.
2416d9e3a6SLisandro Dalcin */
2516d9e3a6SLisandro Dalcin typedef struct {
2616d9e3a6SLisandro Dalcin   HYPRE_Solver       hsolver;
2716d9e3a6SLisandro Dalcin   HYPRE_IJMatrix     ij;
2816d9e3a6SLisandro Dalcin   HYPRE_IJVector     b,x;
2916d9e3a6SLisandro Dalcin 
3016d9e3a6SLisandro Dalcin   PetscErrorCode     (*destroy)(HYPRE_Solver);
3116d9e3a6SLisandro Dalcin   PetscErrorCode     (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
3216d9e3a6SLisandro Dalcin   PetscErrorCode     (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
3316d9e3a6SLisandro Dalcin 
3416d9e3a6SLisandro Dalcin   MPI_Comm           comm_hypre;
3516d9e3a6SLisandro Dalcin   char              *hypre_type;
3616d9e3a6SLisandro Dalcin 
3716d9e3a6SLisandro Dalcin   /* options for Pilut and BoomerAMG*/
3816d9e3a6SLisandro Dalcin   int                maxiter;
3916d9e3a6SLisandro Dalcin   double             tol;
4016d9e3a6SLisandro Dalcin 
4116d9e3a6SLisandro Dalcin   /* options for Pilut */
4216d9e3a6SLisandro Dalcin   int                factorrowsize;
4316d9e3a6SLisandro Dalcin 
4416d9e3a6SLisandro Dalcin   /* options for ParaSails */
4516d9e3a6SLisandro Dalcin   int                nlevels;
4616d9e3a6SLisandro Dalcin   double             threshhold;
4716d9e3a6SLisandro Dalcin   double             filter;
4816d9e3a6SLisandro Dalcin   int                sym;
4916d9e3a6SLisandro Dalcin   double             loadbal;
5016d9e3a6SLisandro Dalcin   int                logging;
5116d9e3a6SLisandro Dalcin   int                ruse;
5216d9e3a6SLisandro Dalcin   int                symt;
5316d9e3a6SLisandro Dalcin 
5416d9e3a6SLisandro Dalcin   /* options for Euclid */
5516d9e3a6SLisandro Dalcin   PetscTruth         bjilu;
5616d9e3a6SLisandro Dalcin   int                levels;
5716d9e3a6SLisandro Dalcin 
5816d9e3a6SLisandro Dalcin   /* options for Euclid and BoomerAMG */
5916d9e3a6SLisandro Dalcin   PetscTruth         printstatistics;
6016d9e3a6SLisandro Dalcin 
6116d9e3a6SLisandro Dalcin   /* options for BoomerAMG */
6216d9e3a6SLisandro Dalcin   int                cycletype;
6316d9e3a6SLisandro Dalcin   int                maxlevels;
6416d9e3a6SLisandro Dalcin   double             strongthreshold;
6516d9e3a6SLisandro Dalcin   double             maxrowsum;
660f1074feSSatish Balay   int                gridsweeps[3];
6716d9e3a6SLisandro Dalcin   int                coarsentype;
6816d9e3a6SLisandro Dalcin   int                measuretype;
690f1074feSSatish Balay   int                relaxtype[3];
7016d9e3a6SLisandro Dalcin   double             relaxweight;
7116d9e3a6SLisandro Dalcin   double             outerrelaxweight;
7216d9e3a6SLisandro Dalcin   int                relaxorder;
7316d9e3a6SLisandro Dalcin   double             truncfactor;
7416d9e3a6SLisandro Dalcin   PetscTruth         applyrichardson;
750f1074feSSatish Balay   int                pmax;
760f1074feSSatish Balay   int                interptype;
770f1074feSSatish Balay   int                agg_nl;
780f1074feSSatish Balay   int                agg_num_paths;
798f87f92bSBarry Smith   int                nodal_coarsen;
808f87f92bSBarry Smith   PetscTruth         nodal_relax;
818f87f92bSBarry Smith   int                nodal_relax_levels;
8216d9e3a6SLisandro Dalcin } PC_HYPRE;
8316d9e3a6SLisandro Dalcin 
8416d9e3a6SLisandro Dalcin 
8516d9e3a6SLisandro Dalcin #undef __FUNCT__
8616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetUp_HYPRE"
8716d9e3a6SLisandro Dalcin static PetscErrorCode PCSetUp_HYPRE(PC pc)
8816d9e3a6SLisandro Dalcin {
8916d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
9016d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
9116d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
9216d9e3a6SLisandro Dalcin   HYPRE_ParVector    bv,xv;
9316d9e3a6SLisandro Dalcin   PetscInt           bs;
9416d9e3a6SLisandro Dalcin   int                hierr;
9516d9e3a6SLisandro Dalcin 
9616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9716d9e3a6SLisandro Dalcin   if (!jac->hypre_type) {
9802a17cd4SBarry Smith     ierr = PCHYPRESetType(pc,"boomeramg");CHKERRQ(ierr);
9916d9e3a6SLisandro Dalcin   }
1005f5c5b43SBarry Smith 
1015f5c5b43SBarry Smith   if (pc->setupcalled) {
1025f5c5b43SBarry Smith     /* always destroy the old matrix and create a new memory;
1035f5c5b43SBarry Smith        hope this does not churn the memory too much. The problem
1045f5c5b43SBarry Smith        is I do not know if it is possible to put the matrix back to
1055f5c5b43SBarry Smith        its initial state so that we can directly copy the values
1065f5c5b43SBarry Smith        the second time through. */
10716d9e3a6SLisandro Dalcin     ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr);
1085f5c5b43SBarry Smith     jac->ij = 0;
10916d9e3a6SLisandro Dalcin   }
1105f5c5b43SBarry Smith 
11116d9e3a6SLisandro Dalcin   if (!jac->ij) { /* create the matrix the first time through */
11216d9e3a6SLisandro Dalcin     ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr);
11316d9e3a6SLisandro Dalcin   }
11416d9e3a6SLisandro Dalcin   if (!jac->b) { /* create the vectors the first time through */
11516d9e3a6SLisandro Dalcin     Vec x,b;
11616d9e3a6SLisandro Dalcin     ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr);
11716d9e3a6SLisandro Dalcin     ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr);
11816d9e3a6SLisandro Dalcin     ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr);
11916d9e3a6SLisandro Dalcin     ierr = VecDestroy(x);CHKERRQ(ierr);
12016d9e3a6SLisandro Dalcin     ierr = VecDestroy(b);CHKERRQ(ierr);
12116d9e3a6SLisandro Dalcin   }
1225f5c5b43SBarry Smith 
12316d9e3a6SLisandro Dalcin   /* special case for BoomerAMG */
12416d9e3a6SLisandro Dalcin   if (jac->setup == HYPRE_BoomerAMGSetup) {
12516d9e3a6SLisandro Dalcin     ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr);
12616d9e3a6SLisandro Dalcin     if (bs > 1) {
12716d9e3a6SLisandro Dalcin       ierr = HYPRE_BoomerAMGSetNumFunctions(jac->hsolver,bs);CHKERRQ(ierr);
12816d9e3a6SLisandro Dalcin     }
12916d9e3a6SLisandro Dalcin   };
13016d9e3a6SLisandro Dalcin   ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr);
13116d9e3a6SLisandro Dalcin   ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr);
13216d9e3a6SLisandro Dalcin   ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&bv);CHKERRQ(ierr);
13316d9e3a6SLisandro Dalcin   ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&xv);CHKERRQ(ierr);
13416d9e3a6SLisandro Dalcin   hierr = (*jac->setup)(jac->hsolver,hmat,bv,xv);
135e32f2f54SBarry Smith   if (hierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE setup, error code %d",hierr);
13616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
13716d9e3a6SLisandro Dalcin }
13816d9e3a6SLisandro Dalcin 
13916d9e3a6SLisandro Dalcin /*
14016d9e3a6SLisandro Dalcin     Replaces the address where the HYPRE vector points to its data with the address of
14116d9e3a6SLisandro Dalcin   PETSc's data. Saves the old address so it can be reset when we are finished with it.
14216d9e3a6SLisandro Dalcin   Allows use to get the data into a HYPRE vector without the cost of memcopies
14316d9e3a6SLisandro Dalcin */
14416d9e3a6SLisandro Dalcin #define HYPREReplacePointer(b,newvalue,savedvalue) {\
14516d9e3a6SLisandro Dalcin    hypre_ParVector *par_vector   = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));\
14616d9e3a6SLisandro Dalcin    hypre_Vector    *local_vector = hypre_ParVectorLocalVector(par_vector);\
14716d9e3a6SLisandro Dalcin    savedvalue         = local_vector->data;\
14816d9e3a6SLisandro Dalcin    local_vector->data = newvalue;}
14916d9e3a6SLisandro Dalcin 
15016d9e3a6SLisandro Dalcin #undef __FUNCT__
15116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApply_HYPRE"
15216d9e3a6SLisandro Dalcin static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
15316d9e3a6SLisandro Dalcin {
15416d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
15516d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
15616d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
15716d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
15816d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
15916d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
16016d9e3a6SLisandro Dalcin   int                hierr;
16116d9e3a6SLisandro Dalcin 
16216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
16316d9e3a6SLisandro Dalcin   if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);}
16416d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
16516d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
16616d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
16716d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
16816d9e3a6SLisandro Dalcin 
16916d9e3a6SLisandro Dalcin   ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr);
17016d9e3a6SLisandro Dalcin   ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr);
17116d9e3a6SLisandro Dalcin   ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr);
17216d9e3a6SLisandro Dalcin   hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv);
1730f1074feSSatish Balay 
174e32f2f54SBarry Smith   /*if (hierr && (hierr != HYPRE_ERROR_CONV || jac->solve != HYPRE_BoomerAMGSolve))SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
1750f1074feSSatish Balay    */
1760f1074feSSatish Balay  /* error code of HYPRE_ERROR_CONV means convergence not achieved - if
1770f1074feSSatish Balay     the tolerance is set to 0.0 (the default), a convergence error will
1780f1074feSSatish Balay     not occur (so we may not want to overide the conv. error here?*/
17965e19b50SBarry Smith  if (hierr && hierr != HYPRE_ERROR_CONV) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
18016d9e3a6SLisandro Dalcin  if (hierr) hypre__global_error = 0;
18116d9e3a6SLisandro Dalcin 
1820f1074feSSatish Balay 
18316d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
18416d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
18516d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
18616d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
18716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
18816d9e3a6SLisandro Dalcin }
18916d9e3a6SLisandro Dalcin 
19016d9e3a6SLisandro Dalcin #undef __FUNCT__
19116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCDestroy_HYPRE"
19216d9e3a6SLisandro Dalcin static PetscErrorCode PCDestroy_HYPRE(PC pc)
19316d9e3a6SLisandro Dalcin {
19416d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
19516d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
19616d9e3a6SLisandro Dalcin 
19716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
19816d9e3a6SLisandro Dalcin   if (jac->ij) { ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); }
19916d9e3a6SLisandro Dalcin   if (jac->b)  { ierr = HYPRE_IJVectorDestroy(jac->b);CHKERRQ(ierr);  }
20016d9e3a6SLisandro Dalcin   if (jac->x)  { ierr = HYPRE_IJVectorDestroy(jac->x);CHKERRQ(ierr);  }
20116d9e3a6SLisandro Dalcin   if (jac->destroy) { ierr = (*jac->destroy)(jac->hsolver);CHKERRQ(ierr); }
202503cfb0cSBarry Smith   ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr);
20316d9e3a6SLisandro Dalcin   if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);}
20416d9e3a6SLisandro Dalcin   ierr = PetscFree(jac);CHKERRQ(ierr);
20516d9e3a6SLisandro Dalcin 
20616d9e3a6SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr);
20716d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);CHKERRQ(ierr);
20816d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);CHKERRQ(ierr);
20916d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
21016d9e3a6SLisandro Dalcin }
21116d9e3a6SLisandro Dalcin 
21216d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
21316d9e3a6SLisandro Dalcin #undef __FUNCT__
21416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut"
21516d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc)
21616d9e3a6SLisandro Dalcin {
21716d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
21816d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
21916d9e3a6SLisandro Dalcin   PetscTruth     flag;
22016d9e3a6SLisandro Dalcin 
22116d9e3a6SLisandro Dalcin   PetscFunctionBegin;
22216d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr);
22316d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr);
22416d9e3a6SLisandro Dalcin   if (flag) {
22516d9e3a6SLisandro Dalcin     ierr = HYPRE_ParCSRPilutSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
22616d9e3a6SLisandro Dalcin   }
22716d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr);
22816d9e3a6SLisandro Dalcin   if (flag) {
22916d9e3a6SLisandro Dalcin     ierr = HYPRE_ParCSRPilutSetDropTolerance(jac->hsolver,jac->tol);CHKERRQ(ierr);
23016d9e3a6SLisandro Dalcin   }
23116d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr);
23216d9e3a6SLisandro Dalcin   if (flag) {
23316d9e3a6SLisandro Dalcin     ierr = HYPRE_ParCSRPilutSetFactorRowSize(jac->hsolver,jac->factorrowsize);CHKERRQ(ierr);
23416d9e3a6SLisandro Dalcin   }
23516d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
23616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
23716d9e3a6SLisandro Dalcin }
23816d9e3a6SLisandro Dalcin 
23916d9e3a6SLisandro Dalcin #undef __FUNCT__
24016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Pilut"
24116d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
24216d9e3a6SLisandro Dalcin {
24316d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
24416d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
24516d9e3a6SLisandro Dalcin   PetscTruth     iascii;
24616d9e3a6SLisandro Dalcin 
24716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
248*2692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
24916d9e3a6SLisandro Dalcin   if (iascii) {
25016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");CHKERRQ(ierr);
25116d9e3a6SLisandro Dalcin     if (jac->maxiter != PETSC_DEFAULT) {
25216d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr);
25316d9e3a6SLisandro Dalcin     } else {
25416d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr);
25516d9e3a6SLisandro Dalcin     }
25616d9e3a6SLisandro Dalcin     if (jac->tol != PETSC_DEFAULT) {
25716d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr);
25816d9e3a6SLisandro Dalcin     } else {
25916d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr);
26016d9e3a6SLisandro Dalcin     }
26116d9e3a6SLisandro Dalcin     if (jac->factorrowsize != PETSC_DEFAULT) {
26216d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr);
26316d9e3a6SLisandro Dalcin     } else {
26416d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");CHKERRQ(ierr);
26516d9e3a6SLisandro Dalcin     }
26616d9e3a6SLisandro Dalcin   }
26716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
26816d9e3a6SLisandro Dalcin }
26916d9e3a6SLisandro Dalcin 
27016d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
27116d9e3a6SLisandro Dalcin #undef __FUNCT__
27216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid"
27316d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc)
27416d9e3a6SLisandro Dalcin {
27516d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
27616d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
27716d9e3a6SLisandro Dalcin   PetscTruth     flag;
278390e7148SBarry Smith   char           *args[8],levels[16];
279390e7148SBarry Smith   PetscInt       cnt = 0;
28016d9e3a6SLisandro Dalcin 
28116d9e3a6SLisandro Dalcin   PetscFunctionBegin;
28216d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr);
28316d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr);
28416d9e3a6SLisandro Dalcin   if (flag) {
28565e19b50SBarry Smith     if (jac->levels < 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels);
28616d9e3a6SLisandro Dalcin     sprintf(levels,"%d",jac->levels);
287390e7148SBarry Smith     args[cnt++] = (char*)"-level"; args[cnt++] = levels;
28816d9e3a6SLisandro Dalcin   }
28916d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);CHKERRQ(ierr);
29016d9e3a6SLisandro Dalcin   if (jac->bjilu) {
291390e7148SBarry Smith     args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1";
29216d9e3a6SLisandro Dalcin   }
29316d9e3a6SLisandro Dalcin 
29416d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);CHKERRQ(ierr);
29516d9e3a6SLisandro Dalcin   if (jac->printstatistics) {
296390e7148SBarry Smith     args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1";
297390e7148SBarry Smith     args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1";
29816d9e3a6SLisandro Dalcin   }
29916d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
300390e7148SBarry Smith   if (cnt) {
301390e7148SBarry Smith     ierr = HYPRE_EuclidSetParams(jac->hsolver,cnt,args);CHKERRQ(ierr);
302390e7148SBarry Smith   }
30316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
30416d9e3a6SLisandro Dalcin }
30516d9e3a6SLisandro Dalcin 
30616d9e3a6SLisandro Dalcin #undef __FUNCT__
30716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Euclid"
30816d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
30916d9e3a6SLisandro Dalcin {
31016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
31116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
31216d9e3a6SLisandro Dalcin   PetscTruth     iascii;
31316d9e3a6SLisandro Dalcin 
31416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
315*2692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
31616d9e3a6SLisandro Dalcin   if (iascii) {
31716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n");CHKERRQ(ierr);
31816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr);
31916d9e3a6SLisandro Dalcin     if (jac->bjilu) {
32016d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr);
32116d9e3a6SLisandro Dalcin     }
32216d9e3a6SLisandro Dalcin   }
32316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
32416d9e3a6SLisandro Dalcin }
32516d9e3a6SLisandro Dalcin 
32616d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
32716d9e3a6SLisandro Dalcin 
32816d9e3a6SLisandro Dalcin #undef __FUNCT__
32916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG"
33016d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
33116d9e3a6SLisandro Dalcin {
33216d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
33316d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
33416d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
33516d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
33616d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
33716d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
33816d9e3a6SLisandro Dalcin   int                hierr;
33916d9e3a6SLisandro Dalcin 
34016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
34116d9e3a6SLisandro Dalcin   ierr = VecSet(x,0.0);CHKERRQ(ierr);
34216d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
34316d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
34416d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
34516d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
34616d9e3a6SLisandro Dalcin 
34716d9e3a6SLisandro Dalcin   ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr);
34816d9e3a6SLisandro Dalcin   ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr);
34916d9e3a6SLisandro Dalcin   ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr);
35016d9e3a6SLisandro Dalcin 
35116d9e3a6SLisandro Dalcin   hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv);
35216d9e3a6SLisandro Dalcin   /* error code of 1 in BoomerAMG merely means convergence not achieved */
353e32f2f54SBarry Smith   if (hierr && (hierr != 1)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
35416d9e3a6SLisandro Dalcin   if (hierr) hypre__global_error = 0;
35516d9e3a6SLisandro Dalcin 
35616d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
35716d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
35816d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
35916d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
36016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
36116d9e3a6SLisandro Dalcin }
36216d9e3a6SLisandro Dalcin 
36316d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
3640f1074feSSatish Balay static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"};
36516d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
3660f1074feSSatish Balay static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","","SOR/Jacobi","backward-SOR/Jacobi","","symmetric-SOR/Jacobi",
3670f1074feSSatish Balay                                                   "","","Gaussian-elimination"};
3680f1074feSSatish Balay static const char *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i",
3690f1074feSSatish Balay                                                   "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"};
37016d9e3a6SLisandro Dalcin #undef __FUNCT__
37116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG"
37216d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
37316d9e3a6SLisandro Dalcin {
37416d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
37516d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
37616d9e3a6SLisandro Dalcin   int            n,indx;
37716d9e3a6SLisandro Dalcin   PetscTruth     flg, tmp_truth;
37816d9e3a6SLisandro Dalcin   double         tmpdbl, twodbl[2];
37916d9e3a6SLisandro Dalcin 
38016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
38116d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr);
38216d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr);
38316d9e3a6SLisandro Dalcin   if (flg) {
38416d9e3a6SLisandro Dalcin     jac->cycletype = indx;
38516d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr);
38616d9e3a6SLisandro Dalcin   }
38716d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr);
38816d9e3a6SLisandro Dalcin   if (flg) {
38965e19b50SBarry Smith     if (jac->maxlevels < 2) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
39016d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr);
39116d9e3a6SLisandro Dalcin   }
39216d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr);
39316d9e3a6SLisandro Dalcin   if (flg) {
39465e19b50SBarry Smith     if (jac->maxiter < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
39516d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
39616d9e3a6SLisandro Dalcin   }
3970f1074feSSatish 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);
39816d9e3a6SLisandro Dalcin   if (flg) {
39965e19b50SBarry 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);
40016d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
40116d9e3a6SLisandro Dalcin   }
40216d9e3a6SLisandro Dalcin 
4030f1074feSSatish Balay   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr);
40416d9e3a6SLisandro Dalcin   if (flg) {
40565e19b50SBarry 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);
40616d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr);
40716d9e3a6SLisandro Dalcin   }
40816d9e3a6SLisandro Dalcin 
4090f1074feSSatish 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);
4100f1074feSSatish Balay   if (flg) {
41165e19b50SBarry 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);
4120f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);CHKERRQ(ierr);
4130f1074feSSatish Balay   }
4140f1074feSSatish Balay 
4150f1074feSSatish Balay  ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);CHKERRQ(ierr);
4160f1074feSSatish Balay   if (flg) {
41765e19b50SBarry 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);
4180f1074feSSatish Balay 
4190f1074feSSatish Balay      ierr = HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl);CHKERRQ(ierr);
4200f1074feSSatish Balay   }
4210f1074feSSatish Balay 
4220f1074feSSatish Balay 
4230f1074feSSatish 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);
4240f1074feSSatish Balay   if (flg) {
42565e19b50SBarry 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);
4260f1074feSSatish Balay 
4270f1074feSSatish Balay      ierr = HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);CHKERRQ(ierr);
4280f1074feSSatish Balay   }
4290f1074feSSatish Balay 
4300f1074feSSatish Balay 
43116d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr);
43216d9e3a6SLisandro Dalcin   if (flg) {
43365e19b50SBarry 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);
43416d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr);
43516d9e3a6SLisandro Dalcin   }
43616d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr);
43716d9e3a6SLisandro Dalcin   if (flg) {
43865e19b50SBarry 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);
43965e19b50SBarry 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);
44016d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr);
44116d9e3a6SLisandro Dalcin   }
44216d9e3a6SLisandro Dalcin 
44316d9e3a6SLisandro Dalcin   /* Grid sweeps */
4440f1074feSSatish 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);
44516d9e3a6SLisandro Dalcin   if (flg) {
44616d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver,indx);CHKERRQ(ierr);
44716d9e3a6SLisandro Dalcin     /* modify the jac structure so we can view the updated options with PC_View */
44816d9e3a6SLisandro Dalcin     jac->gridsweeps[0] = indx;
4490f1074feSSatish Balay     jac->gridsweeps[1] = indx;
4500f1074feSSatish Balay     /*defaults coarse to 1 */
4510f1074feSSatish Balay     jac->gridsweeps[2] = 1;
45216d9e3a6SLisandro Dalcin   }
4530f1074feSSatish Balay 
4540f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx ,&flg);CHKERRQ(ierr);
45516d9e3a6SLisandro Dalcin   if (flg) {
45616d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 1);CHKERRQ(ierr);
4570f1074feSSatish Balay     jac->gridsweeps[0] = indx;
45816d9e3a6SLisandro Dalcin   }
45916d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr);
46016d9e3a6SLisandro Dalcin   if (flg) {
46116d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 2);CHKERRQ(ierr);
4620f1074feSSatish Balay     jac->gridsweeps[1] = indx;
46316d9e3a6SLisandro Dalcin   }
4640f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr);
46516d9e3a6SLisandro Dalcin   if (flg) {
46616d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 3);CHKERRQ(ierr);
4670f1074feSSatish Balay     jac->gridsweeps[2] = indx;
46816d9e3a6SLisandro Dalcin   }
46916d9e3a6SLisandro Dalcin 
47016d9e3a6SLisandro Dalcin   /* Relax type */
4710f1074feSSatish 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);
47216d9e3a6SLisandro Dalcin   if (flg) {
4730f1074feSSatish Balay     jac->relaxtype[0] = jac->relaxtype[1]  = indx;
47416d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, indx);CHKERRQ(ierr);
4750f1074feSSatish Balay     /* by default, coarse type set to 9 */
4760f1074feSSatish Balay     jac->relaxtype[2] = 9;
4770f1074feSSatish Balay 
47816d9e3a6SLisandro Dalcin   }
4790f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
48016d9e3a6SLisandro Dalcin   if (flg) {
48116d9e3a6SLisandro Dalcin     jac->relaxtype[0] = indx;
48216d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 1);CHKERRQ(ierr);
48316d9e3a6SLisandro Dalcin   }
4840f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
48516d9e3a6SLisandro Dalcin   if (flg) {
4860f1074feSSatish Balay     jac->relaxtype[1] = indx;
48716d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 2);CHKERRQ(ierr);
48816d9e3a6SLisandro Dalcin   }
48916d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr);
49016d9e3a6SLisandro Dalcin   if (flg) {
4910f1074feSSatish Balay     jac->relaxtype[2] = indx;
49216d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 3);CHKERRQ(ierr);
49316d9e3a6SLisandro Dalcin   }
49416d9e3a6SLisandro Dalcin 
49516d9e3a6SLisandro Dalcin   /* Relaxation Weight */
49616d9e3a6SLisandro 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);
49716d9e3a6SLisandro Dalcin   if (flg) {
49816d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetRelaxWt(jac->hsolver,tmpdbl);CHKERRQ(ierr);
49916d9e3a6SLisandro Dalcin     jac->relaxweight = tmpdbl;
50016d9e3a6SLisandro Dalcin   }
50116d9e3a6SLisandro Dalcin 
50216d9e3a6SLisandro Dalcin   n=2;
50316d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
50416d9e3a6SLisandro 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);
50516d9e3a6SLisandro Dalcin   if (flg) {
50616d9e3a6SLisandro Dalcin     if (n == 2) {
50716d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
50816d9e3a6SLisandro Dalcin       ierr = HYPRE_BoomerAMGSetLevelRelaxWt(jac->hsolver,twodbl[0],indx);CHKERRQ(ierr);
50916d9e3a6SLisandro Dalcin     } else {
51065e19b50SBarry 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);
51116d9e3a6SLisandro Dalcin     }
51216d9e3a6SLisandro Dalcin   }
51316d9e3a6SLisandro Dalcin 
51416d9e3a6SLisandro Dalcin   /* Outer relaxation Weight */
51516d9e3a6SLisandro 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);
51616d9e3a6SLisandro Dalcin   if (flg) {
51716d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetOuterWt( jac->hsolver, tmpdbl);CHKERRQ(ierr);
51816d9e3a6SLisandro Dalcin     jac->outerrelaxweight = tmpdbl;
51916d9e3a6SLisandro Dalcin   }
52016d9e3a6SLisandro Dalcin 
52116d9e3a6SLisandro Dalcin   n=2;
52216d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
52316d9e3a6SLisandro 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);
52416d9e3a6SLisandro Dalcin   if (flg) {
52516d9e3a6SLisandro Dalcin     if (n == 2) {
52616d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
52716d9e3a6SLisandro Dalcin       ierr = HYPRE_BoomerAMGSetLevelOuterWt( jac->hsolver, twodbl[0], indx);CHKERRQ(ierr);
52816d9e3a6SLisandro Dalcin     } else {
52965e19b50SBarry 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);
53016d9e3a6SLisandro Dalcin     }
53116d9e3a6SLisandro Dalcin   }
53216d9e3a6SLisandro Dalcin 
53316d9e3a6SLisandro Dalcin   /* the Relax Order */
53416d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
53516d9e3a6SLisandro Dalcin 
53616d9e3a6SLisandro Dalcin   if (flg) {
53716d9e3a6SLisandro Dalcin     jac->relaxorder = 0;
53816d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr);
53916d9e3a6SLisandro Dalcin   }
54016d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr);
54116d9e3a6SLisandro Dalcin   if (flg) {
54216d9e3a6SLisandro Dalcin     jac->measuretype = indx;
54316d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr);
54416d9e3a6SLisandro Dalcin   }
5450f1074feSSatish Balay   /* update list length 3/07 */
5460f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,11,HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr);
54716d9e3a6SLisandro Dalcin   if (flg) {
54816d9e3a6SLisandro Dalcin     jac->coarsentype = indx;
54916d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr);
55016d9e3a6SLisandro Dalcin   }
5510f1074feSSatish Balay 
5520f1074feSSatish Balay   /* new 3/07 */
5530f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,14,HYPREBoomerAMGInterpType[0],&indx,&flg);CHKERRQ(ierr);
5540f1074feSSatish Balay   if (flg) {
5550f1074feSSatish Balay     jac->interptype = indx;
5560f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);CHKERRQ(ierr);
5570f1074feSSatish Balay   }
5580f1074feSSatish Balay 
55990d69ab7SBarry Smith   flg  = PETSC_FALSE;
56090d69ab7SBarry Smith   ierr = PetscOptionsTruth("-pc_hypre_boomeramg_print_statistics","Print statistics","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
56116d9e3a6SLisandro Dalcin   if (flg) {
56216d9e3a6SLisandro Dalcin     int level=3;
56316d9e3a6SLisandro Dalcin     jac->printstatistics = PETSC_TRUE;
56416d9e3a6SLisandro Dalcin     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
56516d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetPrintLevel(jac->hsolver,level);CHKERRQ(ierr);
5662ae77aedSBarry Smith   }
5672ae77aedSBarry Smith 
56890d69ab7SBarry Smith   flg  = PETSC_FALSE;
56990d69ab7SBarry Smith   ierr = PetscOptionsTruth("-pc_hypre_boomeramg_print_debug","Print debug information","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
5702ae77aedSBarry Smith   if (flg) {
5712ae77aedSBarry Smith     int level=3;
5722ae77aedSBarry Smith     jac->printstatistics = PETSC_TRUE;
5732ae77aedSBarry Smith     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
57416d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetDebugFlag(jac->hsolver,level);CHKERRQ(ierr);
57516d9e3a6SLisandro Dalcin   }
5768f87f92bSBarry Smith 
5778f87f92bSBarry Smith   ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5788f87f92bSBarry Smith   if (flg && tmp_truth) {
5798f87f92bSBarry Smith     jac->nodal_coarsen = 1;
5808f87f92bSBarry Smith     ierr = HYPRE_BoomerAMGSetNodal(jac->hsolver,1);CHKERRQ(ierr);
5818f87f92bSBarry Smith   }
5828f87f92bSBarry Smith 
5838f87f92bSBarry Smith   ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5848f87f92bSBarry Smith   if (flg && tmp_truth) {
5858f87f92bSBarry Smith     PetscInt tmp_int;
5868f87f92bSBarry Smith     ierr = PetscOptionsInt( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);CHKERRQ(ierr);
5878f87f92bSBarry Smith     if (flg) jac->nodal_relax_levels = tmp_int;
5888f87f92bSBarry Smith     ierr = HYPRE_BoomerAMGSetSmoothType(jac->hsolver,6);CHKERRQ(ierr);
5898f87f92bSBarry Smith     ierr = HYPRE_BoomerAMGSetDomainType(jac->hsolver,1);CHKERRQ(ierr);
5908f87f92bSBarry Smith     ierr = HYPRE_BoomerAMGSetOverlap(jac->hsolver,0);CHKERRQ(ierr);
5918f87f92bSBarry Smith     ierr = HYPRE_BoomerAMGSetSmoothNumLevels(jac->hsolver,jac->nodal_relax_levels);CHKERRQ(ierr);
5928f87f92bSBarry Smith   }
5938f87f92bSBarry Smith 
59416d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
59516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
59616d9e3a6SLisandro Dalcin }
59716d9e3a6SLisandro Dalcin 
59816d9e3a6SLisandro Dalcin #undef __FUNCT__
59916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG"
6007319c654SBarry Smith static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
60116d9e3a6SLisandro Dalcin {
60216d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
60316d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
6044d0a8057SBarry Smith   int            oits;
60516d9e3a6SLisandro Dalcin 
60616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
60716d9e3a6SLisandro Dalcin   ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,its*jac->maxiter);CHKERRQ(ierr);
6084d0a8057SBarry Smith   ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,rtol);CHKERRQ(ierr);
60916d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_TRUE;
61016d9e3a6SLisandro Dalcin   ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr);
61116d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_FALSE;
6124d0a8057SBarry Smith   ierr = HYPRE_BoomerAMGGetNumIterations(jac->hsolver,&oits);CHKERRQ(ierr);
6134d0a8057SBarry Smith   *outits = oits;
6144d0a8057SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
6154d0a8057SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
61616d9e3a6SLisandro Dalcin   ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
61716d9e3a6SLisandro Dalcin   ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
61816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
61916d9e3a6SLisandro Dalcin }
62016d9e3a6SLisandro Dalcin 
62116d9e3a6SLisandro Dalcin 
62216d9e3a6SLisandro Dalcin #undef __FUNCT__
62316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_BoomerAMG"
62416d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
62516d9e3a6SLisandro Dalcin {
62616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
62716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
62816d9e3a6SLisandro Dalcin   PetscTruth     iascii;
62916d9e3a6SLisandro Dalcin 
63016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
631*2692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
63216d9e3a6SLisandro Dalcin   if (iascii) {
63316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr);
63416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr);
63516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr);
63616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr);
63716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr);
63816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr);
6390f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation truncation factor %G\n",jac->truncfactor);CHKERRQ(ierr);
6400f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);CHKERRQ(ierr);
6410f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);CHKERRQ(ierr);
6420f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);CHKERRQ(ierr);
6430f1074feSSatish Balay 
64416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr);
64516d9e3a6SLisandro Dalcin 
6460f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps down         %d\n",jac->gridsweeps[0]);CHKERRQ(ierr);
6470f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps up           %d\n",jac->gridsweeps[1]);CHKERRQ(ierr);
6480f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on coarse    %d\n",jac->gridsweeps[2]);CHKERRQ(ierr);
64916d9e3a6SLisandro Dalcin 
6500f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr);
6510f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr);
6520f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr);
65316d9e3a6SLisandro Dalcin 
65416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax weight  (all)      %G\n",jac->relaxweight);CHKERRQ(ierr);
65516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr);
65616d9e3a6SLisandro Dalcin 
65716d9e3a6SLisandro Dalcin     if (jac->relaxorder) {
65816d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr);
65916d9e3a6SLisandro Dalcin     } else {
66016d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr);
66116d9e3a6SLisandro Dalcin     }
66216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr);
66316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr);
6640f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation type  %s\n",HYPREBoomerAMGInterpType[jac->interptype]);CHKERRQ(ierr);
6658f87f92bSBarry Smith     if (jac->nodal_coarsen) {
6668f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal())\n");CHKERRQ(ierr);
6678f87f92bSBarry Smith     }
6688f87f92bSBarry Smith     if (jac->nodal_relax) {
6698f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);CHKERRQ(ierr);
6708f87f92bSBarry Smith     }
67116d9e3a6SLisandro Dalcin   }
67216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
67316d9e3a6SLisandro Dalcin }
67416d9e3a6SLisandro Dalcin 
67516d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
67616d9e3a6SLisandro Dalcin #undef __FUNCT__
67716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails"
67816d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc)
67916d9e3a6SLisandro Dalcin {
68016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
68116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
68216d9e3a6SLisandro Dalcin   int            indx;
68316d9e3a6SLisandro Dalcin   PetscTruth     flag;
68416d9e3a6SLisandro Dalcin   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};
68516d9e3a6SLisandro Dalcin 
68616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
68716d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr);
68816d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr);
68916d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr);
69016d9e3a6SLisandro Dalcin   if (flag) {
69116d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr);
69216d9e3a6SLisandro Dalcin   }
69316d9e3a6SLisandro Dalcin 
69416d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr);
69516d9e3a6SLisandro Dalcin   if (flag) {
69616d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr);
69716d9e3a6SLisandro Dalcin   }
69816d9e3a6SLisandro Dalcin 
69916d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr);
70016d9e3a6SLisandro Dalcin   if (flag) {
70116d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr);
70216d9e3a6SLisandro Dalcin   }
70316d9e3a6SLisandro Dalcin 
70416d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,&flag);CHKERRQ(ierr);
70516d9e3a6SLisandro Dalcin   if (flag) {
70616d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr);
70716d9e3a6SLisandro Dalcin   }
70816d9e3a6SLisandro Dalcin 
70916d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,&flag);CHKERRQ(ierr);
71016d9e3a6SLisandro Dalcin   if (flag) {
71116d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr);
71216d9e3a6SLisandro Dalcin   }
71316d9e3a6SLisandro Dalcin 
71416d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);CHKERRQ(ierr);
71516d9e3a6SLisandro Dalcin   if (flag) {
71616d9e3a6SLisandro Dalcin     jac->symt = indx;
71716d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr);
71816d9e3a6SLisandro Dalcin   }
71916d9e3a6SLisandro Dalcin 
72016d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
72116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
72216d9e3a6SLisandro Dalcin }
72316d9e3a6SLisandro Dalcin 
72416d9e3a6SLisandro Dalcin #undef __FUNCT__
72516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_ParaSails"
72616d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
72716d9e3a6SLisandro Dalcin {
72816d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
72916d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
73016d9e3a6SLisandro Dalcin   PetscTruth     iascii;
73116d9e3a6SLisandro Dalcin   const char     *symt = 0;;
73216d9e3a6SLisandro Dalcin 
73316d9e3a6SLisandro Dalcin   PetscFunctionBegin;
734*2692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
73516d9e3a6SLisandro Dalcin   if (iascii) {
73616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");CHKERRQ(ierr);
73716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr);
73816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr);
73916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr);
74016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr);
74116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscTruths[jac->ruse]);CHKERRQ(ierr);
74216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscTruths[jac->logging]);CHKERRQ(ierr);
74316d9e3a6SLisandro Dalcin     if (!jac->symt) {
74416d9e3a6SLisandro Dalcin       symt = "nonsymmetric matrix and preconditioner";
74516d9e3a6SLisandro Dalcin     } else if (jac->symt == 1) {
74616d9e3a6SLisandro Dalcin       symt = "SPD matrix and preconditioner";
74716d9e3a6SLisandro Dalcin     } else if (jac->symt == 2) {
74816d9e3a6SLisandro Dalcin       symt = "nonsymmetric matrix but SPD preconditioner";
74916d9e3a6SLisandro Dalcin     } else {
75065e19b50SBarry Smith       SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
75116d9e3a6SLisandro Dalcin     }
75216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr);
75316d9e3a6SLisandro Dalcin   }
75416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
75516d9e3a6SLisandro Dalcin }
75616d9e3a6SLisandro Dalcin /* ---------------------------------------------------------------------------------*/
75716d9e3a6SLisandro Dalcin 
75816d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
75916d9e3a6SLisandro Dalcin #undef __FUNCT__
76016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType_HYPRE"
76116d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType_HYPRE(PC pc,const char *name[])
76216d9e3a6SLisandro Dalcin {
76316d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
76416d9e3a6SLisandro Dalcin 
76516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
76616d9e3a6SLisandro Dalcin   *name = jac->hypre_type;
76716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
76816d9e3a6SLisandro Dalcin }
76916d9e3a6SLisandro Dalcin EXTERN_C_END
77016d9e3a6SLisandro Dalcin 
77116d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
77216d9e3a6SLisandro Dalcin #undef __FUNCT__
77316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType_HYPRE"
77416d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType_HYPRE(PC pc,const char name[])
77516d9e3a6SLisandro Dalcin {
77616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
77716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
77816d9e3a6SLisandro Dalcin   PetscTruth     flag;
77916d9e3a6SLisandro Dalcin 
78016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
78116d9e3a6SLisandro Dalcin   if (jac->hypre_type) {
78216d9e3a6SLisandro Dalcin     ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr);
783e7e72b3dSBarry Smith     if (!flag) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
78416d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
78516d9e3a6SLisandro Dalcin   } else {
78616d9e3a6SLisandro Dalcin     ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr);
78716d9e3a6SLisandro Dalcin   }
78816d9e3a6SLisandro Dalcin 
78916d9e3a6SLisandro Dalcin   jac->maxiter            = PETSC_DEFAULT;
79016d9e3a6SLisandro Dalcin   jac->tol                = PETSC_DEFAULT;
79116d9e3a6SLisandro Dalcin   jac->printstatistics    = PetscLogPrintInfo;
79216d9e3a6SLisandro Dalcin 
79316d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr);
79416d9e3a6SLisandro Dalcin   if (flag) {
79516d9e3a6SLisandro Dalcin     ierr                    = HYPRE_ParCSRPilutCreate(jac->comm_hypre,&jac->hsolver);
79616d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
79716d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Pilut;
79816d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParCSRPilutDestroy;
79916d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParCSRPilutSetup;
80016d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParCSRPilutSolve;
80116d9e3a6SLisandro Dalcin     jac->factorrowsize      = PETSC_DEFAULT;
80216d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
80316d9e3a6SLisandro Dalcin   }
80416d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr);
80516d9e3a6SLisandro Dalcin   if (flag) {
80616d9e3a6SLisandro Dalcin     ierr                    = HYPRE_ParaSailsCreate(jac->comm_hypre,&jac->hsolver);
80716d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
80816d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_ParaSails;
80916d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParaSailsDestroy;
81016d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParaSailsSetup;
81116d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParaSailsSolve;
81216d9e3a6SLisandro Dalcin     /* initialize */
81316d9e3a6SLisandro Dalcin     jac->nlevels     = 1;
81416d9e3a6SLisandro Dalcin     jac->threshhold  = .1;
81516d9e3a6SLisandro Dalcin     jac->filter      = .1;
81616d9e3a6SLisandro Dalcin     jac->loadbal     = 0;
81716d9e3a6SLisandro Dalcin     if (PetscLogPrintInfo) {
81816d9e3a6SLisandro Dalcin       jac->logging     = (int) PETSC_TRUE;
81916d9e3a6SLisandro Dalcin     } else {
82016d9e3a6SLisandro Dalcin       jac->logging     = (int) PETSC_FALSE;
82116d9e3a6SLisandro Dalcin     }
82216d9e3a6SLisandro Dalcin     jac->ruse = (int) PETSC_FALSE;
82316d9e3a6SLisandro Dalcin     jac->symt   = 0;
82416d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr);
82516d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr);
82616d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr);
82716d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr);
82816d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr);
82916d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr);
83016d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
83116d9e3a6SLisandro Dalcin   }
83216d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr);
83316d9e3a6SLisandro Dalcin   if (flag) {
83416d9e3a6SLisandro Dalcin     ierr                    = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver);
83516d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
83616d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Euclid;
83716d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_EuclidDestroy;
83816d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_EuclidSetup;
83916d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_EuclidSolve;
84016d9e3a6SLisandro Dalcin     /* initialization */
84116d9e3a6SLisandro Dalcin     jac->bjilu              = PETSC_FALSE;
84216d9e3a6SLisandro Dalcin     jac->levels             = 1;
84316d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
84416d9e3a6SLisandro Dalcin   }
84516d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr);
84616d9e3a6SLisandro Dalcin   if (flag) {
84716d9e3a6SLisandro Dalcin     ierr                     = HYPRE_BoomerAMGCreate(&jac->hsolver);
84816d9e3a6SLisandro Dalcin     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
84916d9e3a6SLisandro Dalcin     pc->ops->view            = PCView_HYPRE_BoomerAMG;
85016d9e3a6SLisandro Dalcin     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
85116d9e3a6SLisandro Dalcin     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
85216d9e3a6SLisandro Dalcin     jac->destroy             = HYPRE_BoomerAMGDestroy;
85316d9e3a6SLisandro Dalcin     jac->setup               = HYPRE_BoomerAMGSetup;
85416d9e3a6SLisandro Dalcin     jac->solve               = HYPRE_BoomerAMGSolve;
85516d9e3a6SLisandro Dalcin     jac->applyrichardson     = PETSC_FALSE;
85616d9e3a6SLisandro Dalcin     /* these defaults match the hypre defaults */
85716d9e3a6SLisandro Dalcin     jac->cycletype        = 1;
85816d9e3a6SLisandro Dalcin     jac->maxlevels        = 25;
85916d9e3a6SLisandro Dalcin     jac->maxiter          = 1;
8608f87f92bSBarry Smith     jac->tol              = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */
86116d9e3a6SLisandro Dalcin     jac->truncfactor      = 0.0;
86216d9e3a6SLisandro Dalcin     jac->strongthreshold  = .25;
86316d9e3a6SLisandro Dalcin     jac->maxrowsum        = .9;
86416d9e3a6SLisandro Dalcin     jac->coarsentype      = 6;
86516d9e3a6SLisandro Dalcin     jac->measuretype      = 0;
8660f1074feSSatish Balay     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = 1;
8678f87f92bSBarry Smith     jac->relaxtype[0]     = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */
8680f1074feSSatish Balay     jac->relaxtype[2]     = 9; /*G.E. */
86916d9e3a6SLisandro Dalcin     jac->relaxweight      = 1.0;
87016d9e3a6SLisandro Dalcin     jac->outerrelaxweight = 1.0;
87116d9e3a6SLisandro Dalcin     jac->relaxorder       = 1;
8720f1074feSSatish Balay     jac->interptype       = 0;
8730f1074feSSatish Balay     jac->agg_nl           = 0;
8740f1074feSSatish Balay     jac->pmax             = 0;
8750f1074feSSatish Balay     jac->truncfactor      = 0.0;
8760f1074feSSatish Balay     jac->agg_num_paths    = 1;
8778f87f92bSBarry Smith 
8788f87f92bSBarry Smith     jac->nodal_coarsen    = 0;
8798f87f92bSBarry Smith     jac->nodal_relax      = PETSC_FALSE;
8808f87f92bSBarry Smith     jac->nodal_relax_levels = 1;
88116d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr);
88216d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr);
88316d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
88416d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
88516d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr);
88616d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr);
88716d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr);
88816d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr);
88916d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr);
89016d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr);
8910f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);CHKERRQ(ierr);
8920f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl);
8930f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);CHKERRQ(ierr);
8940f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);CHKERRQ(ierr);
8950f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, jac->relaxtype[0]);  /*defaults coarse to 9*/
8960f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver, jac->gridsweeps[0]); /*defaults coarse to 1 */
89716d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
89816d9e3a6SLisandro Dalcin   }
899503cfb0cSBarry Smith   ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr);
90016d9e3a6SLisandro Dalcin   jac->hypre_type = PETSC_NULL;
90165e19b50SBarry Smith   SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
90216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
90316d9e3a6SLisandro Dalcin }
90416d9e3a6SLisandro Dalcin EXTERN_C_END
90516d9e3a6SLisandro Dalcin 
90616d9e3a6SLisandro Dalcin /*
90716d9e3a6SLisandro Dalcin     It only gets here if the HYPRE type has not been set before the call to
90816d9e3a6SLisandro Dalcin    ...SetFromOptions() which actually is most of the time
90916d9e3a6SLisandro Dalcin */
91016d9e3a6SLisandro Dalcin #undef __FUNCT__
91116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE"
91216d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE(PC pc)
91316d9e3a6SLisandro Dalcin {
91416d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
91516d9e3a6SLisandro Dalcin   int            indx;
91616d9e3a6SLisandro Dalcin   const char     *type[] = {"pilut","parasails","boomeramg","euclid"};
91716d9e3a6SLisandro Dalcin   PetscTruth     flg;
91816d9e3a6SLisandro Dalcin 
91916d9e3a6SLisandro Dalcin   PetscFunctionBegin;
92016d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr);
92102a17cd4SBarry Smith   ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"boomeramg",&indx,&flg);CHKERRQ(ierr);
92216d9e3a6SLisandro Dalcin   if (flg) {
92316d9e3a6SLisandro Dalcin     ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr);
92402a17cd4SBarry Smith   } else {
92502a17cd4SBarry Smith     ierr = PCHYPRESetType_HYPRE(pc,"boomeramg");CHKERRQ(ierr);
92616d9e3a6SLisandro Dalcin   }
92716d9e3a6SLisandro Dalcin   if (pc->ops->setfromoptions) {
92816d9e3a6SLisandro Dalcin     ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr);
92916d9e3a6SLisandro Dalcin   }
93016d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
93116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
93216d9e3a6SLisandro Dalcin }
93316d9e3a6SLisandro Dalcin 
93416d9e3a6SLisandro Dalcin #undef __FUNCT__
93516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType"
93616d9e3a6SLisandro Dalcin /*@C
93716d9e3a6SLisandro Dalcin      PCHYPRESetType - Sets which hypre preconditioner you wish to use
93816d9e3a6SLisandro Dalcin 
93916d9e3a6SLisandro Dalcin    Input Parameters:
94016d9e3a6SLisandro Dalcin +     pc - the preconditioner context
94116d9e3a6SLisandro Dalcin -     name - either  pilut, parasails, boomeramg, euclid
94216d9e3a6SLisandro Dalcin 
94316d9e3a6SLisandro Dalcin    Options Database Keys:
94416d9e3a6SLisandro Dalcin    -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
94516d9e3a6SLisandro Dalcin 
94616d9e3a6SLisandro Dalcin    Level: intermediate
94716d9e3a6SLisandro Dalcin 
94816d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
94916d9e3a6SLisandro Dalcin            PCHYPRE
95016d9e3a6SLisandro Dalcin 
95116d9e3a6SLisandro Dalcin @*/
95216d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC pc,const char name[])
95316d9e3a6SLisandro Dalcin {
95416d9e3a6SLisandro Dalcin   PetscErrorCode ierr,(*f)(PC,const char[]);
95516d9e3a6SLisandro Dalcin 
95616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9570700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
95816d9e3a6SLisandro Dalcin   PetscValidCharPointer(name,2);
95916d9e3a6SLisandro Dalcin   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPRESetType_C",(void (**)(void))&f);CHKERRQ(ierr);
96016d9e3a6SLisandro Dalcin   if (f) {
96116d9e3a6SLisandro Dalcin     ierr = (*f)(pc,name);CHKERRQ(ierr);
96216d9e3a6SLisandro Dalcin   }
96316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
96416d9e3a6SLisandro Dalcin }
96516d9e3a6SLisandro Dalcin 
96616d9e3a6SLisandro Dalcin #undef __FUNCT__
96716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType"
96816d9e3a6SLisandro Dalcin /*@C
96916d9e3a6SLisandro Dalcin      PCHYPREGetType - Gets which hypre preconditioner you are using
97016d9e3a6SLisandro Dalcin 
97116d9e3a6SLisandro Dalcin    Input Parameter:
97216d9e3a6SLisandro Dalcin .     pc - the preconditioner context
97316d9e3a6SLisandro Dalcin 
97416d9e3a6SLisandro Dalcin    Output Parameter:
97516d9e3a6SLisandro Dalcin .     name - either  pilut, parasails, boomeramg, euclid
97616d9e3a6SLisandro Dalcin 
97716d9e3a6SLisandro Dalcin    Level: intermediate
97816d9e3a6SLisandro Dalcin 
97916d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC,
98016d9e3a6SLisandro Dalcin            PCHYPRE
98116d9e3a6SLisandro Dalcin 
98216d9e3a6SLisandro Dalcin @*/
98316d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC pc,const char *name[])
98416d9e3a6SLisandro Dalcin {
98516d9e3a6SLisandro Dalcin   PetscErrorCode ierr,(*f)(PC,const char*[]);
98616d9e3a6SLisandro Dalcin 
98716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9880700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
98916d9e3a6SLisandro Dalcin   PetscValidPointer(name,2);
99016d9e3a6SLisandro Dalcin   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPREGetType_C",(void (**)(void))&f);CHKERRQ(ierr);
99116d9e3a6SLisandro Dalcin   if (f) {
99216d9e3a6SLisandro Dalcin     ierr = (*f)(pc,name);CHKERRQ(ierr);
99316d9e3a6SLisandro Dalcin   }
99416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
99516d9e3a6SLisandro Dalcin }
99616d9e3a6SLisandro Dalcin 
99716d9e3a6SLisandro Dalcin /*MC
99816d9e3a6SLisandro Dalcin      PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre
99916d9e3a6SLisandro Dalcin 
100016d9e3a6SLisandro Dalcin    Options Database Keys:
100116d9e3a6SLisandro Dalcin +   -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
100216d9e3a6SLisandro Dalcin -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
100316d9e3a6SLisandro Dalcin           preconditioner
100416d9e3a6SLisandro Dalcin 
100516d9e3a6SLisandro Dalcin    Level: intermediate
100616d9e3a6SLisandro Dalcin 
100716d9e3a6SLisandro Dalcin    Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()),
100816d9e3a6SLisandro Dalcin           the many hypre options can ONLY be set via the options database (e.g. the command line
100916d9e3a6SLisandro Dalcin           or with PetscOptionsSetValue(), there are no functions to set them)
101016d9e3a6SLisandro Dalcin 
101116d9e3a6SLisandro Dalcin           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations
10120f1074feSSatish Balay           (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if
10130f1074feSSatish Balay           -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner
10140f1074feSSatish Balay           (-pc_hypre_boomeramg_rtol should be set to 0.0 - the default - to strictly use a fixed number of
10158f87f92bSBarry Smith           iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations
10160f1074feSSatish Balay           and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10
10170f1074feSSatish Balay           then AT MOST twenty V-cycles of boomeramg will be called.
101816d9e3a6SLisandro Dalcin 
10190f1074feSSatish Balay            Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation
10200f1074feSSatish Balay            (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry.
10210f1074feSSatish Balay            Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi.
102216d9e3a6SLisandro Dalcin           If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly
102316d9e3a6SLisandro Dalcin           and use -ksp_max_it to control the number of V-cycles.
102416d9e3a6SLisandro Dalcin           (see the PETSc FAQ.html at the PETSc website under the Documentation tab).
102516d9e3a6SLisandro Dalcin 
102616d9e3a6SLisandro Dalcin 	  2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option
102716d9e3a6SLisandro Dalcin 	  -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L.
102816d9e3a6SLisandro Dalcin 
10299e5bc791SBarry Smith           See PCPFMG for access to the hypre Struct PFMG solver
10309e5bc791SBarry Smith 
103116d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
10329e5bc791SBarry Smith            PCHYPRESetType(), PCPFMG
103316d9e3a6SLisandro Dalcin 
103416d9e3a6SLisandro Dalcin M*/
103516d9e3a6SLisandro Dalcin 
103616d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
103716d9e3a6SLisandro Dalcin #undef __FUNCT__
103816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCCreate_HYPRE"
103916d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_HYPRE(PC pc)
104016d9e3a6SLisandro Dalcin {
104116d9e3a6SLisandro Dalcin   PC_HYPRE       *jac;
104216d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
104316d9e3a6SLisandro Dalcin 
104416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
104538f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_HYPRE,&jac);CHKERRQ(ierr);
104616d9e3a6SLisandro Dalcin   pc->data                 = jac;
104716d9e3a6SLisandro Dalcin   pc->ops->destroy         = PCDestroy_HYPRE;
104816d9e3a6SLisandro Dalcin   pc->ops->setfromoptions  = PCSetFromOptions_HYPRE;
104916d9e3a6SLisandro Dalcin   pc->ops->setup           = PCSetUp_HYPRE;
105016d9e3a6SLisandro Dalcin   pc->ops->apply           = PCApply_HYPRE;
105116d9e3a6SLisandro Dalcin   jac->comm_hypre          = MPI_COMM_NULL;
105216d9e3a6SLisandro Dalcin   jac->hypre_type          = PETSC_NULL;
105316d9e3a6SLisandro Dalcin   /* duplicate communicator for hypre */
10547adad957SLisandro Dalcin   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(jac->comm_hypre));CHKERRQ(ierr);
105516d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr);
105616d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr);
105716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
105816d9e3a6SLisandro Dalcin }
105916d9e3a6SLisandro Dalcin EXTERN_C_END
1060ebc551c0SBarry Smith 
1061f91d8e95SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------*/
1062f91d8e95SBarry Smith 
1063b862ddfaSBarry Smith /* working with a HYPRE_StructMatrix and need access to its data */
1064ee691612SBarry Smith #include "../src/dm/da/utils/mhyp.h"
1065b862ddfaSBarry Smith /* this include is needed ONLY to allow access to the private data inside the Mat object specific to hypre */
106668326731SBarry Smith #include "private/matimpl.h"
1067ebc551c0SBarry Smith 
1068ebc551c0SBarry Smith typedef struct {
106968326731SBarry Smith   MPI_Comm            hcomm;       /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */
1070f91d8e95SBarry Smith   HYPRE_StructSolver  hsolver;
10719e5bc791SBarry Smith 
10729e5bc791SBarry Smith   /* keep copy of PFMG options used so may view them */
10739e5bc791SBarry Smith   int                 its;
10749e5bc791SBarry Smith   double              tol;
10759e5bc791SBarry Smith   int                 relax_type;
10769e5bc791SBarry Smith   int                 rap_type;
10779e5bc791SBarry Smith   int                 num_pre_relax,num_post_relax;
10783b46a515SGlenn Hammond   int                 max_levels;
1079ebc551c0SBarry Smith } PC_PFMG;
1080ebc551c0SBarry Smith 
1081ebc551c0SBarry Smith #undef __FUNCT__
1082ebc551c0SBarry Smith #define __FUNCT__ "PCDestroy_PFMG"
1083ebc551c0SBarry Smith PetscErrorCode PCDestroy_PFMG(PC pc)
1084ebc551c0SBarry Smith {
1085ebc551c0SBarry Smith   PetscErrorCode ierr;
1086f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1087ebc551c0SBarry Smith 
1088ebc551c0SBarry Smith   PetscFunctionBegin;
1089f91d8e95SBarry Smith   if (ex->hsolver) {ierr = HYPRE_StructPFMGDestroy(ex->hsolver);CHKERRQ(ierr);}
1090f91d8e95SBarry Smith   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1091f91d8e95SBarry Smith   ierr = PetscFree(ex);CHKERRQ(ierr);
1092ebc551c0SBarry Smith   PetscFunctionReturn(0);
1093ebc551c0SBarry Smith }
1094ebc551c0SBarry Smith 
10959e5bc791SBarry Smith static const char *PFMGRelaxType[]   = {"Jacobi","Weighted-Jacobi","symmetric-Red/Black-Gauss-Seidel","Red/Black-Gauss-Seidel"};
10969e5bc791SBarry Smith static const char *PFMGRAPType[]   = {"Galerkin","non-Galerkin"};
10979e5bc791SBarry Smith 
1098ebc551c0SBarry Smith #undef __FUNCT__
1099ebc551c0SBarry Smith #define __FUNCT__ "PCView_PFMG"
1100ebc551c0SBarry Smith PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer)
1101ebc551c0SBarry Smith {
1102ebc551c0SBarry Smith   PetscErrorCode ierr;
1103ebc551c0SBarry Smith   PetscTruth     iascii;
1104f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1105ebc551c0SBarry Smith 
1106ebc551c0SBarry Smith   PetscFunctionBegin;
1107*2692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
11089e5bc791SBarry Smith   if (iascii) {
11099e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG preconditioning\n");CHKERRQ(ierr);
11109e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
11119e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
11129e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
11139e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: RAP type %s\n",PFMGRAPType[ex->rap_type]);CHKERRQ(ierr);
11149e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
11153b46a515SGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max levels %d\n",ex->max_levels);CHKERRQ(ierr);
11169e5bc791SBarry Smith   }
1117ebc551c0SBarry Smith   PetscFunctionReturn(0);
1118ebc551c0SBarry Smith }
1119ebc551c0SBarry Smith 
11209e5bc791SBarry Smith 
1121ebc551c0SBarry Smith #undef __FUNCT__
1122ebc551c0SBarry Smith #define __FUNCT__ "PCSetFromOptions_PFMG"
1123ebc551c0SBarry Smith PetscErrorCode PCSetFromOptions_PFMG(PC pc)
1124ebc551c0SBarry Smith {
1125ebc551c0SBarry Smith   PetscErrorCode ierr;
1126f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1127886061d4SBarry Smith   PetscTruth     flg = PETSC_FALSE;
1128ebc551c0SBarry Smith 
1129ebc551c0SBarry Smith   PetscFunctionBegin;
1130ebc551c0SBarry Smith   ierr = PetscOptionsHead("PFMG options");CHKERRQ(ierr);
11319e5bc791SBarry Smith   ierr = PetscOptionsTruth("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
113268326731SBarry Smith   if (flg) {
1133a0324ebeSBarry Smith     int level=3;
113468326731SBarry Smith     ierr = HYPRE_StructPFMGSetPrintLevel(ex->hsolver,level);CHKERRQ(ierr);
113568326731SBarry Smith   }
11369e5bc791SBarry 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);
11379e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetMaxIter(ex->hsolver,ex->its);CHKERRQ(ierr);
11389e5bc791SBarry 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);
11399e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetNumPreRelax(ex->hsolver,ex->num_pre_relax);CHKERRQ(ierr);
11409e5bc791SBarry 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);
11419e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetNumPostRelax(ex->hsolver,ex->num_post_relax);CHKERRQ(ierr);
11429e5bc791SBarry Smith 
11433b46a515SGlenn Hammond   ierr = PetscOptionsInt("-pc_pfmg_max_levels","Max Levels for MG hierarchy","HYPRE_StructPFMGSetMaxLevels",ex->max_levels,&ex->max_levels,PETSC_NULL);CHKERRQ(ierr);
11443b46a515SGlenn Hammond   ierr = HYPRE_StructPFMGSetMaxLevels(ex->hsolver,ex->max_levels);CHKERRQ(ierr);
11453b46a515SGlenn Hammond 
11469e5bc791SBarry Smith   ierr = PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr);
11479e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetTol(ex->hsolver,ex->tol);CHKERRQ(ierr);
11489e5bc791SBarry 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);
11499e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetRelaxType(ex->hsolver, ex->relax_type);CHKERRQ(ierr);
11509e5bc791SBarry Smith   ierr = PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,2,PFMGRAPType[ex->rap_type],&ex->rap_type,PETSC_NULL);CHKERRQ(ierr);
11519e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetRAPType(ex->hsolver, ex->rap_type);CHKERRQ(ierr);
1152ebc551c0SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1153ebc551c0SBarry Smith   PetscFunctionReturn(0);
1154ebc551c0SBarry Smith }
1155ebc551c0SBarry Smith 
1156f91d8e95SBarry Smith #undef __FUNCT__
1157f91d8e95SBarry Smith #define __FUNCT__ "PCApply_PFMG"
1158f91d8e95SBarry Smith PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y)
1159f91d8e95SBarry Smith {
1160f91d8e95SBarry Smith   PetscErrorCode  ierr;
1161f91d8e95SBarry Smith   PC_PFMG         *ex = (PC_PFMG*) pc->data;
1162f91d8e95SBarry Smith   PetscScalar     *xx,*yy;
1163f91d8e95SBarry Smith   int             ilower[3],iupper[3];
116468326731SBarry Smith   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
1165f91d8e95SBarry Smith 
1166f91d8e95SBarry Smith   PetscFunctionBegin;
116768326731SBarry Smith   ierr = DAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1168f91d8e95SBarry Smith   iupper[0] += ilower[0] - 1;
1169f91d8e95SBarry Smith   iupper[1] += ilower[1] - 1;
1170f91d8e95SBarry Smith   iupper[2] += ilower[2] - 1;
1171f91d8e95SBarry Smith 
1172f91d8e95SBarry Smith   /* copy x values over to hypre */
1173c7ef6ce2SGlenn Hammond   ierr = HYPRE_StructVectorSetConstantValues(mx->hb,0.0);CHKERRQ(ierr);
1174f91d8e95SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
117568326731SBarry Smith   ierr = HYPRE_StructVectorSetBoxValues(mx->hb,ilower,iupper,xx);CHKERRQ(ierr);
1176f91d8e95SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
117768326731SBarry Smith   ierr = HYPRE_StructVectorAssemble(mx->hb);CHKERRQ(ierr);
1178f91d8e95SBarry Smith 
117968326731SBarry Smith   ierr = HYPRE_StructPFMGSolve(ex->hsolver,mx->hmat,mx->hb,mx->hx);CHKERRQ(ierr);
1180f91d8e95SBarry Smith 
1181f91d8e95SBarry Smith   /* copy solution values back to PETSc */
1182f91d8e95SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1183a0324ebeSBarry Smith   ierr = HYPRE_StructVectorGetBoxValues(mx->hx,ilower,iupper,yy);CHKERRQ(ierr);
1184f91d8e95SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1185f91d8e95SBarry Smith   PetscFunctionReturn(0);
1186f91d8e95SBarry Smith }
1187f91d8e95SBarry Smith 
11889e5bc791SBarry Smith #undef __FUNCT__
11899e5bc791SBarry Smith #define __FUNCT__ "PCApplyRichardson_PFMG"
11907319c654SBarry Smith static PetscErrorCode PCApplyRichardson_PFMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
11919e5bc791SBarry Smith {
11929e5bc791SBarry Smith   PC_PFMG        *jac = (PC_PFMG*)pc->data;
11939e5bc791SBarry Smith   PetscErrorCode ierr;
11949e5bc791SBarry Smith   int            oits;
11959e5bc791SBarry Smith 
11969e5bc791SBarry Smith   PetscFunctionBegin;
11979e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetMaxIter(jac->hsolver,its*jac->its);CHKERRQ(ierr);
11989e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetTol(jac->hsolver,rtol);CHKERRQ(ierr);
11999e5bc791SBarry Smith 
12009e5bc791SBarry Smith   ierr = PCApply_PFMG(pc,b,y);CHKERRQ(ierr);
12019e5bc791SBarry Smith   ierr = HYPRE_StructPFMGGetNumIterations(jac->hsolver,&oits);CHKERRQ(ierr);
12029e5bc791SBarry Smith   *outits = oits;
12039e5bc791SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
12049e5bc791SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
12059e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
12069e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetMaxIter(jac->hsolver,jac->its);CHKERRQ(ierr);
12079e5bc791SBarry Smith   PetscFunctionReturn(0);
12089e5bc791SBarry Smith }
12099e5bc791SBarry Smith 
12109e5bc791SBarry Smith 
12113a32d3dbSGlenn Hammond #undef __FUNCT__
12123a32d3dbSGlenn Hammond #define __FUNCT__ "PCSetUp_PFMG"
12133a32d3dbSGlenn Hammond PetscErrorCode PCSetUp_PFMG(PC pc)
12143a32d3dbSGlenn Hammond {
12153a32d3dbSGlenn Hammond   PetscErrorCode  ierr;
12163a32d3dbSGlenn Hammond   PC_PFMG         *ex = (PC_PFMG*) pc->data;
12173a32d3dbSGlenn Hammond   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
12183a32d3dbSGlenn Hammond   PetscTruth      flg;
12193a32d3dbSGlenn Hammond 
12203a32d3dbSGlenn Hammond   PetscFunctionBegin;
12213a32d3dbSGlenn Hammond   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATHYPRESTRUCT,&flg);CHKERRQ(ierr);
1222e7e72b3dSBarry Smith   if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESTRUCT with this preconditioner");
12233a32d3dbSGlenn Hammond 
12243a32d3dbSGlenn Hammond   /* create the hypre solver object and set its information */
12253a32d3dbSGlenn Hammond   if (ex->hsolver) {
12263a32d3dbSGlenn Hammond     ierr = HYPRE_StructPFMGDestroy(ex->hsolver);CHKERRQ(ierr);
12273a32d3dbSGlenn Hammond   }
12283a32d3dbSGlenn Hammond   ierr = HYPRE_StructPFMGCreate(ex->hcomm,&ex->hsolver);CHKERRQ(ierr);
12293a32d3dbSGlenn Hammond   ierr = PCSetFromOptions_PFMG(pc);CHKERRQ(ierr);
12303a32d3dbSGlenn Hammond   ierr = HYPRE_StructPFMGSetup(ex->hsolver,mx->hmat,mx->hb,mx->hx);CHKERRQ(ierr);
12313a32d3dbSGlenn Hammond   ierr = HYPRE_StructPFMGSetZeroGuess(ex->hsolver);CHKERRQ(ierr);
12323a32d3dbSGlenn Hammond 
12333a32d3dbSGlenn Hammond   PetscFunctionReturn(0);
12343a32d3dbSGlenn Hammond }
12353a32d3dbSGlenn Hammond 
1236ebc551c0SBarry Smith 
1237ebc551c0SBarry Smith /*MC
1238ebc551c0SBarry Smith      PCPFMG - the hypre PFMG multigrid solver
1239ebc551c0SBarry Smith 
1240ebc551c0SBarry Smith    Level: advanced
1241ebc551c0SBarry Smith 
12429e5bc791SBarry Smith    Options Database:
12439e5bc791SBarry Smith + -pc_pfmg_its <its> number of iterations of PFMG to use as preconditioner
12449e5bc791SBarry Smith . -pc_pfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
12459e5bc791SBarry Smith . -pc_pfmg_num_post_relax <steps> number of smoothing steps after coarse grid
12469e5bc791SBarry Smith . -pc_pfmg_tol <tol> tolerance of PFMG
12479e5bc791SBarry 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
12489e5bc791SBarry Smith - -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin
1249f91d8e95SBarry Smith 
12509e5bc791SBarry Smith    Notes:  This is for CELL-centered descretizations
12519e5bc791SBarry Smith 
12528e395302SJed Brown            This must be used with the MATHYPRESTRUCT matrix type.
12539e5bc791SBarry Smith            This is less general than in hypre, it supports only one block per process defined by a PETSc DA.
12549e5bc791SBarry Smith 
12559e5bc791SBarry Smith .seealso:  PCMG, MATHYPRESTRUCT
1256ebc551c0SBarry Smith M*/
1257ebc551c0SBarry Smith 
1258ebc551c0SBarry Smith EXTERN_C_BEGIN
1259ebc551c0SBarry Smith #undef __FUNCT__
1260ebc551c0SBarry Smith #define __FUNCT__ "PCCreate_PFMG"
1261ebc551c0SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_PFMG(PC pc)
1262ebc551c0SBarry Smith {
1263ebc551c0SBarry Smith   PetscErrorCode ierr;
1264ebc551c0SBarry Smith   PC_PFMG        *ex;
1265ebc551c0SBarry Smith 
1266ebc551c0SBarry Smith   PetscFunctionBegin;
1267ebc551c0SBarry Smith   ierr = PetscNew(PC_PFMG,&ex);CHKERRQ(ierr);\
126868326731SBarry Smith   pc->data = ex;
1269ebc551c0SBarry Smith 
12709e5bc791SBarry Smith   ex->its            = 1;
12719e5bc791SBarry Smith   ex->tol            = 1.e-8;
12729e5bc791SBarry Smith   ex->relax_type     = 1;
12739e5bc791SBarry Smith   ex->rap_type       = 0;
12749e5bc791SBarry Smith   ex->num_pre_relax  = 1;
12759e5bc791SBarry Smith   ex->num_post_relax = 1;
12763b46a515SGlenn Hammond   ex->max_levels     = 0;
12779e5bc791SBarry Smith 
1278ebc551c0SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_PFMG;
1279ebc551c0SBarry Smith   pc->ops->view            = PCView_PFMG;
1280ebc551c0SBarry Smith   pc->ops->destroy         = PCDestroy_PFMG;
1281f91d8e95SBarry Smith   pc->ops->apply           = PCApply_PFMG;
12829e5bc791SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_PFMG;
128368326731SBarry Smith   pc->ops->setup           = PCSetUp_PFMG;
1284f91d8e95SBarry Smith   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr);
128568326731SBarry Smith   ierr = HYPRE_StructPFMGCreate(ex->hcomm,&ex->hsolver);CHKERRQ(ierr);
1286ebc551c0SBarry Smith   PetscFunctionReturn(0);
1287ebc551c0SBarry Smith }
1288ebc551c0SBarry Smith EXTERN_C_END
1289d851a50bSGlenn Hammond 
1290d851a50bSGlenn Hammond /* we know we are working with a HYPRE_SStructMatrix */
1291d851a50bSGlenn Hammond typedef struct {
1292d851a50bSGlenn Hammond   MPI_Comm            hcomm;       /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */
1293d851a50bSGlenn Hammond   HYPRE_SStructSolver  ss_solver;
1294d851a50bSGlenn Hammond 
1295d851a50bSGlenn Hammond   /* keep copy of SYSPFMG options used so may view them */
1296d851a50bSGlenn Hammond   int                 its;
1297d851a50bSGlenn Hammond   double              tol;
1298d851a50bSGlenn Hammond   int                 relax_type;
1299d851a50bSGlenn Hammond   int                 num_pre_relax,num_post_relax;
1300d851a50bSGlenn Hammond } PC_SysPFMG;
1301d851a50bSGlenn Hammond 
1302d851a50bSGlenn Hammond #undef __FUNCT__
1303d851a50bSGlenn Hammond #define __FUNCT__ "PCDestroy_SysPFMG"
1304d851a50bSGlenn Hammond PetscErrorCode PCDestroy_SysPFMG(PC pc)
1305d851a50bSGlenn Hammond {
1306d851a50bSGlenn Hammond   PetscErrorCode ierr;
1307d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1308d851a50bSGlenn Hammond 
1309d851a50bSGlenn Hammond   PetscFunctionBegin;
1310d851a50bSGlenn Hammond   if (ex->ss_solver) {ierr = HYPRE_SStructSysPFMGDestroy(ex->ss_solver);CHKERRQ(ierr);}
1311d851a50bSGlenn Hammond   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1312d851a50bSGlenn Hammond   ierr = PetscFree(ex);CHKERRQ(ierr);
1313d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1314d851a50bSGlenn Hammond }
1315d851a50bSGlenn Hammond 
1316d851a50bSGlenn Hammond static const char *SysPFMGRelaxType[]   = {"Weighted-Jacobi","Red/Black-Gauss-Seidel"};
1317d851a50bSGlenn Hammond 
1318d851a50bSGlenn Hammond #undef __FUNCT__
1319d851a50bSGlenn Hammond #define __FUNCT__ "PCView_SysPFMG"
1320d851a50bSGlenn Hammond PetscErrorCode PCView_SysPFMG(PC pc,PetscViewer viewer)
1321d851a50bSGlenn Hammond {
1322d851a50bSGlenn Hammond   PetscErrorCode ierr;
1323d851a50bSGlenn Hammond   PetscTruth     iascii;
1324d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1325d851a50bSGlenn Hammond 
1326d851a50bSGlenn Hammond   PetscFunctionBegin;
1327*2692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1328d851a50bSGlenn Hammond   if (iascii) {
1329d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG preconditioning\n");CHKERRQ(ierr);
1330d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
1331d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
1332d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
1333d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
1334d851a50bSGlenn Hammond   }
1335d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1336d851a50bSGlenn Hammond }
1337d851a50bSGlenn Hammond 
1338d851a50bSGlenn Hammond 
1339d851a50bSGlenn Hammond #undef __FUNCT__
1340d851a50bSGlenn Hammond #define __FUNCT__ "PCSetFromOptions_SysPFMG"
1341d851a50bSGlenn Hammond PetscErrorCode PCSetFromOptions_SysPFMG(PC pc)
1342d851a50bSGlenn Hammond {
1343d851a50bSGlenn Hammond   PetscErrorCode ierr;
1344d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1345d851a50bSGlenn Hammond   PetscTruth     flg = PETSC_FALSE;
1346d851a50bSGlenn Hammond 
1347d851a50bSGlenn Hammond   PetscFunctionBegin;
1348d851a50bSGlenn Hammond   ierr = PetscOptionsHead("SysPFMG options");CHKERRQ(ierr);
1349d851a50bSGlenn Hammond   ierr = PetscOptionsTruth("-pc_syspfmg_print_statistics","Print statistics","HYPRE_SStructSysPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1350d851a50bSGlenn Hammond   if (flg) {
1351d851a50bSGlenn Hammond     int level=3;
1352d851a50bSGlenn Hammond     ierr = HYPRE_SStructSysPFMGSetPrintLevel(ex->ss_solver,level);CHKERRQ(ierr);
1353d851a50bSGlenn Hammond   }
1354d851a50bSGlenn 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);
1355d851a50bSGlenn Hammond   ierr = HYPRE_SStructSysPFMGSetMaxIter(ex->ss_solver,ex->its);CHKERRQ(ierr);
1356d851a50bSGlenn 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);
1357d851a50bSGlenn Hammond   ierr = HYPRE_SStructSysPFMGSetNumPreRelax(ex->ss_solver,ex->num_pre_relax);CHKERRQ(ierr);
1358d851a50bSGlenn 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);
1359d851a50bSGlenn Hammond   ierr = HYPRE_SStructSysPFMGSetNumPostRelax(ex->ss_solver,ex->num_post_relax);CHKERRQ(ierr);
1360d851a50bSGlenn Hammond 
1361d851a50bSGlenn Hammond   ierr = PetscOptionsReal("-pc_syspfmg_tol","Tolerance of SysPFMG","HYPRE_SStructSysPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr);
1362d851a50bSGlenn Hammond   ierr = HYPRE_SStructSysPFMGSetTol(ex->ss_solver,ex->tol);CHKERRQ(ierr);
1363d851a50bSGlenn 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);
1364d851a50bSGlenn Hammond   ierr = HYPRE_SStructSysPFMGSetRelaxType(ex->ss_solver, ex->relax_type);CHKERRQ(ierr);
1365d851a50bSGlenn Hammond   ierr = PetscOptionsTail();CHKERRQ(ierr);
1366d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1367d851a50bSGlenn Hammond }
1368d851a50bSGlenn Hammond 
1369d851a50bSGlenn Hammond #undef __FUNCT__
1370d851a50bSGlenn Hammond #define __FUNCT__ "PCApply_SysPFMG"
1371d851a50bSGlenn Hammond PetscErrorCode PCApply_SysPFMG(PC pc,Vec x,Vec y)
1372d851a50bSGlenn Hammond {
1373d851a50bSGlenn Hammond   PetscErrorCode    ierr;
1374d851a50bSGlenn Hammond   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
1375d851a50bSGlenn Hammond   PetscScalar      *xx,*yy;
1376d851a50bSGlenn Hammond   int               ilower[3],iupper[3];
1377d851a50bSGlenn Hammond   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct *)(pc->pmat->data);
1378d851a50bSGlenn Hammond   int               ordering= mx->dofs_order;
1379d851a50bSGlenn Hammond   int               nvars= mx->nvars;
1380d851a50bSGlenn Hammond   int               part= 0;
1381d851a50bSGlenn Hammond   int               size;
1382d851a50bSGlenn Hammond   int               i;
1383d851a50bSGlenn Hammond 
1384d851a50bSGlenn Hammond   PetscFunctionBegin;
1385d851a50bSGlenn Hammond   ierr = DAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1386d851a50bSGlenn Hammond   iupper[0] += ilower[0] - 1;
1387d851a50bSGlenn Hammond   iupper[1] += ilower[1] - 1;
1388d851a50bSGlenn Hammond   iupper[2] += ilower[2] - 1;
1389d851a50bSGlenn Hammond 
1390d851a50bSGlenn Hammond   size= 1;
1391d851a50bSGlenn Hammond   for (i= 0; i< 3; i++) {
1392d851a50bSGlenn Hammond      size*= (iupper[i]-ilower[i]+1);
1393d851a50bSGlenn Hammond   }
1394d851a50bSGlenn Hammond   /* copy x values over to hypre for variable ordering */
1395d851a50bSGlenn Hammond   if (ordering) {
1396d851a50bSGlenn Hammond      ierr = HYPRE_SStructVectorSetConstantValues(mx->ss_b,0.0);CHKERRQ(ierr);
1397d851a50bSGlenn Hammond      ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1398d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
1399d851a50bSGlenn Hammond         ierr = HYPRE_SStructVectorSetBoxValues(mx->ss_b,part,ilower,iupper,i,xx+(size*i));CHKERRQ(ierr);
1400d851a50bSGlenn Hammond      }
1401d851a50bSGlenn Hammond      ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1402d851a50bSGlenn Hammond      ierr = HYPRE_SStructVectorAssemble(mx->ss_b);CHKERRQ(ierr);
1403d851a50bSGlenn Hammond 
1404d851a50bSGlenn Hammond      ierr = HYPRE_SStructMatrixMatvec(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x);CHKERRQ(ierr);
1405d851a50bSGlenn Hammond      ierr = HYPRE_SStructSysPFMGSolve(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x);CHKERRQ(ierr);
1406d851a50bSGlenn Hammond 
1407d851a50bSGlenn Hammond      /* copy solution values back to PETSc */
1408d851a50bSGlenn Hammond      ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1409d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
1410d851a50bSGlenn Hammond         ierr = HYPRE_SStructVectorGetBoxValues(mx->ss_x,part,ilower,iupper,i,yy+(size*i));CHKERRQ(ierr);
1411d851a50bSGlenn Hammond      }
1412d851a50bSGlenn Hammond      ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1413d851a50bSGlenn Hammond   }
1414d851a50bSGlenn Hammond 
1415d851a50bSGlenn Hammond   else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
1416d851a50bSGlenn Hammond      PetscScalar     *z;
1417d851a50bSGlenn Hammond      int              j, k;
1418d851a50bSGlenn Hammond 
1419d851a50bSGlenn Hammond      ierr = PetscMalloc(nvars*size*sizeof(PetscScalar),&z);CHKERRQ(ierr);
1420d851a50bSGlenn Hammond      ierr = HYPRE_SStructVectorSetConstantValues(mx->ss_b,0.0);CHKERRQ(ierr);
1421d851a50bSGlenn Hammond      ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1422d851a50bSGlenn Hammond 
1423d851a50bSGlenn Hammond      /* transform nodal to hypre's variable ordering for sys_pfmg */
1424d851a50bSGlenn Hammond      for (i= 0; i< size; i++) {
1425d851a50bSGlenn Hammond         k= i*nvars;
1426d851a50bSGlenn Hammond         for (j= 0; j< nvars; j++) {
1427d851a50bSGlenn Hammond            z[j*size+i]= xx[k+j];
1428d851a50bSGlenn Hammond         }
1429d851a50bSGlenn Hammond      }
1430d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
1431d851a50bSGlenn Hammond         ierr = HYPRE_SStructVectorSetBoxValues(mx->ss_b,part,ilower,iupper,i,z+(size*i));CHKERRQ(ierr);
1432d851a50bSGlenn Hammond      }
1433d851a50bSGlenn Hammond      ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1434d851a50bSGlenn Hammond 
1435d851a50bSGlenn Hammond      ierr = HYPRE_SStructVectorAssemble(mx->ss_b);CHKERRQ(ierr);
1436d851a50bSGlenn Hammond 
1437d851a50bSGlenn Hammond      ierr = HYPRE_SStructSysPFMGSolve(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x);CHKERRQ(ierr);
1438d851a50bSGlenn Hammond 
1439d851a50bSGlenn Hammond      /* copy solution values back to PETSc */
1440d851a50bSGlenn Hammond      ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1441d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
1442d851a50bSGlenn Hammond         ierr = HYPRE_SStructVectorGetBoxValues(mx->ss_x,part,ilower,iupper,i,z+(size*i));CHKERRQ(ierr);
1443d851a50bSGlenn Hammond      }
1444d851a50bSGlenn Hammond      /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
1445d851a50bSGlenn Hammond      for (i= 0; i< size; i++) {
1446d851a50bSGlenn Hammond         k= i*nvars;
1447d851a50bSGlenn Hammond         for (j= 0; j< nvars; j++) {
1448d851a50bSGlenn Hammond            yy[k+j]= z[j*size+i];
1449d851a50bSGlenn Hammond         }
1450d851a50bSGlenn Hammond      }
1451d851a50bSGlenn Hammond      ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1452d851a50bSGlenn Hammond 
1453d851a50bSGlenn Hammond      ierr = PetscFree(z);CHKERRQ(ierr);
1454d851a50bSGlenn Hammond   }
1455d851a50bSGlenn Hammond 
1456d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1457d851a50bSGlenn Hammond }
1458d851a50bSGlenn Hammond 
1459d851a50bSGlenn Hammond #undef __FUNCT__
1460d851a50bSGlenn Hammond #define __FUNCT__ "PCApplyRichardson_SysPFMG"
1461d851a50bSGlenn Hammond static PetscErrorCode PCApplyRichardson_SysPFMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
1462d851a50bSGlenn Hammond {
1463d851a50bSGlenn Hammond   PC_SysPFMG    *jac = (PC_SysPFMG*)pc->data;
1464d851a50bSGlenn Hammond   PetscErrorCode ierr;
1465d851a50bSGlenn Hammond   int            oits;
1466d851a50bSGlenn Hammond 
1467d851a50bSGlenn Hammond   PetscFunctionBegin;
1468d851a50bSGlenn Hammond 
1469d851a50bSGlenn Hammond   ierr = HYPRE_SStructSysPFMGSetMaxIter(jac->ss_solver,its*jac->its);CHKERRQ(ierr);
1470d851a50bSGlenn Hammond   ierr = HYPRE_SStructSysPFMGSetTol(jac->ss_solver,rtol);CHKERRQ(ierr);
1471d851a50bSGlenn Hammond 
1472d851a50bSGlenn Hammond   ierr = PCApply_SysPFMG(pc,b,y);CHKERRQ(ierr);
1473d851a50bSGlenn Hammond   ierr = HYPRE_SStructSysPFMGGetNumIterations(jac->ss_solver,&oits);CHKERRQ(ierr);
1474d851a50bSGlenn Hammond   *outits = oits;
1475d851a50bSGlenn Hammond   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
1476d851a50bSGlenn Hammond   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
1477d851a50bSGlenn Hammond   ierr = HYPRE_SStructSysPFMGSetTol(jac->ss_solver,jac->tol);CHKERRQ(ierr);
1478d851a50bSGlenn Hammond   ierr = HYPRE_SStructSysPFMGSetMaxIter(jac->ss_solver,jac->its);CHKERRQ(ierr);
1479d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1480d851a50bSGlenn Hammond }
1481d851a50bSGlenn Hammond 
1482d851a50bSGlenn Hammond 
1483d851a50bSGlenn Hammond #undef __FUNCT__
1484d851a50bSGlenn Hammond #define __FUNCT__ "PCSetUp_SysPFMG"
1485d851a50bSGlenn Hammond PetscErrorCode PCSetUp_SysPFMG(PC pc)
1486d851a50bSGlenn Hammond {
1487d851a50bSGlenn Hammond   PetscErrorCode    ierr;
1488d851a50bSGlenn Hammond   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
1489d851a50bSGlenn Hammond   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct *)(pc->pmat->data);
1490d851a50bSGlenn Hammond   PetscTruth        flg;
1491d851a50bSGlenn Hammond 
1492d851a50bSGlenn Hammond   PetscFunctionBegin;
1493d851a50bSGlenn Hammond   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATHYPRESSTRUCT,&flg);CHKERRQ(ierr);
1494e7e72b3dSBarry Smith   if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESSTRUCT with this preconditioner");
1495d851a50bSGlenn Hammond 
1496d851a50bSGlenn Hammond   /* create the hypre sstruct solver object and set its information */
1497d851a50bSGlenn Hammond   if (ex->ss_solver) {
1498d851a50bSGlenn Hammond     ierr = HYPRE_SStructSysPFMGDestroy(ex->ss_solver);CHKERRQ(ierr);
1499d851a50bSGlenn Hammond   }
1500d851a50bSGlenn Hammond   ierr = HYPRE_SStructSysPFMGCreate(ex->hcomm,&ex->ss_solver);CHKERRQ(ierr);
1501d851a50bSGlenn Hammond   ierr = PCSetFromOptions_SysPFMG(pc);CHKERRQ(ierr);
1502d851a50bSGlenn Hammond   ierr = HYPRE_SStructSysPFMGSetZeroGuess(ex->ss_solver);CHKERRQ(ierr);
1503d851a50bSGlenn Hammond   ierr = HYPRE_SStructSysPFMGSetup(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x);CHKERRQ(ierr);
1504d851a50bSGlenn Hammond 
1505d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1506d851a50bSGlenn Hammond }
1507d851a50bSGlenn Hammond 
1508d851a50bSGlenn Hammond 
1509d851a50bSGlenn Hammond /*MC
1510d851a50bSGlenn Hammond      PCSysPFMG - the hypre SysPFMG multigrid solver
1511d851a50bSGlenn Hammond 
1512d851a50bSGlenn Hammond    Level: advanced
1513d851a50bSGlenn Hammond 
1514d851a50bSGlenn Hammond    Options Database:
1515d851a50bSGlenn Hammond + -pc_syspfmg_its <its> number of iterations of SysPFMG to use as preconditioner
1516d851a50bSGlenn Hammond . -pc_syspfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
1517d851a50bSGlenn Hammond . -pc_syspfmg_num_post_relax <steps> number of smoothing steps after coarse grid
1518d851a50bSGlenn Hammond . -pc_syspfmg_tol <tol> tolerance of SysPFMG
1519d851a50bSGlenn Hammond . -pc_syspfmg_relax_type -relaxation type for the up and down cycles, one of Weighted-Jacobi,Red/Black-Gauss-Seidel
1520d851a50bSGlenn Hammond 
1521d851a50bSGlenn Hammond    Notes:  This is for CELL-centered descretizations
1522d851a50bSGlenn Hammond 
1523f6680f47SSatish Balay            This must be used with the MATHYPRESSTRUCT matrix type.
1524d851a50bSGlenn Hammond            This is less general than in hypre, it supports only one part, and one block per process defined by a PETSc DA.
1525d851a50bSGlenn Hammond            Also, only cell-centered variables.
1526d851a50bSGlenn Hammond 
1527d851a50bSGlenn Hammond .seealso:  PCMG, MATHYPRESSTRUCT
1528d851a50bSGlenn Hammond M*/
1529d851a50bSGlenn Hammond 
1530d851a50bSGlenn Hammond EXTERN_C_BEGIN
1531d851a50bSGlenn Hammond #undef __FUNCT__
1532d851a50bSGlenn Hammond #define __FUNCT__ "PCCreate_SysPFMG"
1533d851a50bSGlenn Hammond PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_SysPFMG(PC pc)
1534d851a50bSGlenn Hammond {
1535d851a50bSGlenn Hammond   PetscErrorCode     ierr;
1536d851a50bSGlenn Hammond   PC_SysPFMG        *ex;
1537d851a50bSGlenn Hammond 
1538d851a50bSGlenn Hammond   PetscFunctionBegin;
1539d851a50bSGlenn Hammond   ierr = PetscNew(PC_SysPFMG,&ex);CHKERRQ(ierr);\
1540d851a50bSGlenn Hammond   pc->data = ex;
1541d851a50bSGlenn Hammond 
1542d851a50bSGlenn Hammond   ex->its            = 1;
1543d851a50bSGlenn Hammond   ex->tol            = 1.e-8;
1544d851a50bSGlenn Hammond   ex->relax_type     = 1;
1545d851a50bSGlenn Hammond   ex->num_pre_relax  = 1;
1546d851a50bSGlenn Hammond   ex->num_post_relax = 1;
1547d851a50bSGlenn Hammond 
1548d851a50bSGlenn Hammond   pc->ops->setfromoptions  = PCSetFromOptions_SysPFMG;
1549d851a50bSGlenn Hammond   pc->ops->view            = PCView_SysPFMG;
1550d851a50bSGlenn Hammond   pc->ops->destroy         = PCDestroy_SysPFMG;
1551d851a50bSGlenn Hammond   pc->ops->apply           = PCApply_SysPFMG;
1552d851a50bSGlenn Hammond   pc->ops->applyrichardson = PCApplyRichardson_SysPFMG;
1553d851a50bSGlenn Hammond   pc->ops->setup           = PCSetUp_SysPFMG;
1554d851a50bSGlenn Hammond   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr);
1555d851a50bSGlenn Hammond   ierr = HYPRE_SStructSysPFMGCreate(ex->hcomm,&ex->ss_solver);CHKERRQ(ierr);
1556d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1557d851a50bSGlenn Hammond }
1558d851a50bSGlenn Hammond EXTERN_C_END
1559