xref: /petsc/src/ksp/pc/impls/hypre/hypre.c (revision c7ef6ce2c37faeeecf5bb5fa630bca9d6fdfea17)
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);
13516d9e3a6SLisandro Dalcin   if (hierr) SETERRQ1(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 
1740f1074feSSatish Balay   /*if (hierr && (hierr != HYPRE_ERROR_CONV || jac->solve != HYPRE_BoomerAMGSolve))SETERRQ1(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?*/
1790f1074feSSatish Balay  if (hierr && hierr != HYPRE_ERROR_CONV)
1800f1074feSSatish Balay   {
1810f1074feSSatish Balay      SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
1820f1074feSSatish Balay   }
18316d9e3a6SLisandro Dalcin  if (hierr) hypre__global_error = 0;
18416d9e3a6SLisandro Dalcin 
1850f1074feSSatish Balay 
18616d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
18716d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
18816d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
18916d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
19016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
19116d9e3a6SLisandro Dalcin }
19216d9e3a6SLisandro Dalcin 
19316d9e3a6SLisandro Dalcin #undef __FUNCT__
19416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCDestroy_HYPRE"
19516d9e3a6SLisandro Dalcin static PetscErrorCode PCDestroy_HYPRE(PC pc)
19616d9e3a6SLisandro Dalcin {
19716d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
19816d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
19916d9e3a6SLisandro Dalcin 
20016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
20116d9e3a6SLisandro Dalcin   if (jac->ij) { ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); }
20216d9e3a6SLisandro Dalcin   if (jac->b)  { ierr = HYPRE_IJVectorDestroy(jac->b);CHKERRQ(ierr);  }
20316d9e3a6SLisandro Dalcin   if (jac->x)  { ierr = HYPRE_IJVectorDestroy(jac->x);CHKERRQ(ierr);  }
20416d9e3a6SLisandro Dalcin   if (jac->destroy) { ierr = (*jac->destroy)(jac->hsolver);CHKERRQ(ierr); }
20516d9e3a6SLisandro Dalcin   ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr);
20616d9e3a6SLisandro Dalcin   if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);}
20716d9e3a6SLisandro Dalcin   ierr = PetscFree(jac);CHKERRQ(ierr);
20816d9e3a6SLisandro Dalcin 
20916d9e3a6SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr);
21016d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);CHKERRQ(ierr);
21116d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);CHKERRQ(ierr);
21216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
21316d9e3a6SLisandro Dalcin }
21416d9e3a6SLisandro Dalcin 
21516d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
21616d9e3a6SLisandro Dalcin #undef __FUNCT__
21716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut"
21816d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc)
21916d9e3a6SLisandro Dalcin {
22016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
22116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
22216d9e3a6SLisandro Dalcin   PetscTruth     flag;
22316d9e3a6SLisandro Dalcin 
22416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
22516d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr);
22616d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr);
22716d9e3a6SLisandro Dalcin   if (flag) {
22816d9e3a6SLisandro Dalcin     ierr = HYPRE_ParCSRPilutSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
22916d9e3a6SLisandro Dalcin   }
23016d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr);
23116d9e3a6SLisandro Dalcin   if (flag) {
23216d9e3a6SLisandro Dalcin     ierr = HYPRE_ParCSRPilutSetDropTolerance(jac->hsolver,jac->tol);CHKERRQ(ierr);
23316d9e3a6SLisandro Dalcin   }
23416d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr);
23516d9e3a6SLisandro Dalcin   if (flag) {
23616d9e3a6SLisandro Dalcin     ierr = HYPRE_ParCSRPilutSetFactorRowSize(jac->hsolver,jac->factorrowsize);CHKERRQ(ierr);
23716d9e3a6SLisandro Dalcin   }
23816d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
23916d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
24016d9e3a6SLisandro Dalcin }
24116d9e3a6SLisandro Dalcin 
24216d9e3a6SLisandro Dalcin #undef __FUNCT__
24316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Pilut"
24416d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
24516d9e3a6SLisandro Dalcin {
24616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
24716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
24816d9e3a6SLisandro Dalcin   PetscTruth     iascii;
24916d9e3a6SLisandro Dalcin 
25016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
25116d9e3a6SLisandro Dalcin   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
25216d9e3a6SLisandro Dalcin   if (iascii) {
25316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");CHKERRQ(ierr);
25416d9e3a6SLisandro Dalcin     if (jac->maxiter != PETSC_DEFAULT) {
25516d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr);
25616d9e3a6SLisandro Dalcin     } else {
25716d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr);
25816d9e3a6SLisandro Dalcin     }
25916d9e3a6SLisandro Dalcin     if (jac->tol != PETSC_DEFAULT) {
26016d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr);
26116d9e3a6SLisandro Dalcin     } else {
26216d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr);
26316d9e3a6SLisandro Dalcin     }
26416d9e3a6SLisandro Dalcin     if (jac->factorrowsize != PETSC_DEFAULT) {
26516d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr);
26616d9e3a6SLisandro Dalcin     } else {
26716d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");CHKERRQ(ierr);
26816d9e3a6SLisandro Dalcin     }
26916d9e3a6SLisandro Dalcin   }
27016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
27116d9e3a6SLisandro Dalcin }
27216d9e3a6SLisandro Dalcin 
27316d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
27416d9e3a6SLisandro Dalcin #undef __FUNCT__
27516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid"
27616d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc)
27716d9e3a6SLisandro Dalcin {
27816d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
27916d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
28016d9e3a6SLisandro Dalcin   PetscTruth     flag;
281390e7148SBarry Smith   char           *args[8],levels[16];
282390e7148SBarry Smith   PetscInt       cnt = 0;
28316d9e3a6SLisandro Dalcin 
28416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
28516d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr);
28616d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr);
28716d9e3a6SLisandro Dalcin   if (flag) {
28816d9e3a6SLisandro Dalcin     if (jac->levels < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels);
28916d9e3a6SLisandro Dalcin     sprintf(levels,"%d",jac->levels);
290390e7148SBarry Smith     args[cnt++] = (char*)"-level"; args[cnt++] = levels;
29116d9e3a6SLisandro Dalcin   }
29216d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);CHKERRQ(ierr);
29316d9e3a6SLisandro Dalcin   if (jac->bjilu) {
294390e7148SBarry Smith     args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1";
29516d9e3a6SLisandro Dalcin   }
29616d9e3a6SLisandro Dalcin 
29716d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);CHKERRQ(ierr);
29816d9e3a6SLisandro Dalcin   if (jac->printstatistics) {
299390e7148SBarry Smith     args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1";
300390e7148SBarry Smith     args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1";
30116d9e3a6SLisandro Dalcin   }
30216d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
303390e7148SBarry Smith   if (cnt) {
304390e7148SBarry Smith     ierr = HYPRE_EuclidSetParams(jac->hsolver,cnt,args);CHKERRQ(ierr);
305390e7148SBarry Smith   }
30616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
30716d9e3a6SLisandro Dalcin }
30816d9e3a6SLisandro Dalcin 
30916d9e3a6SLisandro Dalcin #undef __FUNCT__
31016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Euclid"
31116d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
31216d9e3a6SLisandro Dalcin {
31316d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
31416d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
31516d9e3a6SLisandro Dalcin   PetscTruth     iascii;
31616d9e3a6SLisandro Dalcin 
31716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
31816d9e3a6SLisandro Dalcin   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
31916d9e3a6SLisandro Dalcin   if (iascii) {
32016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n");CHKERRQ(ierr);
32116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr);
32216d9e3a6SLisandro Dalcin     if (jac->bjilu) {
32316d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr);
32416d9e3a6SLisandro Dalcin     }
32516d9e3a6SLisandro Dalcin   }
32616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
32716d9e3a6SLisandro Dalcin }
32816d9e3a6SLisandro Dalcin 
32916d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
33016d9e3a6SLisandro Dalcin 
33116d9e3a6SLisandro Dalcin #undef __FUNCT__
33216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG"
33316d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
33416d9e3a6SLisandro Dalcin {
33516d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
33616d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
33716d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
33816d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
33916d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
34016d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
34116d9e3a6SLisandro Dalcin   int                hierr;
34216d9e3a6SLisandro Dalcin 
34316d9e3a6SLisandro Dalcin   PetscFunctionBegin;
34416d9e3a6SLisandro Dalcin   ierr = VecSet(x,0.0);CHKERRQ(ierr);
34516d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
34616d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
34716d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
34816d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
34916d9e3a6SLisandro Dalcin 
35016d9e3a6SLisandro Dalcin   ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr);
35116d9e3a6SLisandro Dalcin   ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr);
35216d9e3a6SLisandro Dalcin   ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr);
35316d9e3a6SLisandro Dalcin 
35416d9e3a6SLisandro Dalcin   hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv);
35516d9e3a6SLisandro Dalcin   /* error code of 1 in BoomerAMG merely means convergence not achieved */
35616d9e3a6SLisandro Dalcin   if (hierr && (hierr != 1)) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
35716d9e3a6SLisandro Dalcin   if (hierr) hypre__global_error = 0;
35816d9e3a6SLisandro Dalcin 
35916d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
36016d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
36116d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
36216d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
36316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
36416d9e3a6SLisandro Dalcin }
36516d9e3a6SLisandro Dalcin 
36616d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
3670f1074feSSatish Balay static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"};
36816d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
3690f1074feSSatish Balay static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","","SOR/Jacobi","backward-SOR/Jacobi","","symmetric-SOR/Jacobi",
3700f1074feSSatish Balay                                                   "","","Gaussian-elimination"};
3710f1074feSSatish Balay static const char *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i",
3720f1074feSSatish Balay                                                   "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"};
37316d9e3a6SLisandro Dalcin #undef __FUNCT__
37416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG"
37516d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
37616d9e3a6SLisandro Dalcin {
37716d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
37816d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
37916d9e3a6SLisandro Dalcin   int            n,indx;
38016d9e3a6SLisandro Dalcin   PetscTruth     flg, tmp_truth;
38116d9e3a6SLisandro Dalcin   double         tmpdbl, twodbl[2];
38216d9e3a6SLisandro Dalcin 
38316d9e3a6SLisandro Dalcin   PetscFunctionBegin;
38416d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr);
38516d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr);
38616d9e3a6SLisandro Dalcin   if (flg) {
38716d9e3a6SLisandro Dalcin     jac->cycletype = indx;
38816d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr);
38916d9e3a6SLisandro Dalcin   }
39016d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr);
39116d9e3a6SLisandro Dalcin   if (flg) {
39216d9e3a6SLisandro Dalcin     if (jac->maxlevels < 2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
39316d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr);
39416d9e3a6SLisandro Dalcin   }
39516d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr);
39616d9e3a6SLisandro Dalcin   if (flg) {
39716d9e3a6SLisandro Dalcin     if (jac->maxiter < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
39816d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
39916d9e3a6SLisandro Dalcin   }
4000f1074feSSatish 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);
40116d9e3a6SLisandro Dalcin   if (flg) {
4020f1074feSSatish Balay     if (jac->tol < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be greater than or equal to zero",jac->tol);
40316d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
40416d9e3a6SLisandro Dalcin   }
40516d9e3a6SLisandro Dalcin 
4060f1074feSSatish Balay   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr);
40716d9e3a6SLisandro Dalcin   if (flg) {
40816d9e3a6SLisandro Dalcin     if (jac->truncfactor < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor);
40916d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr);
41016d9e3a6SLisandro Dalcin   }
41116d9e3a6SLisandro Dalcin 
4120f1074feSSatish 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);
4130f1074feSSatish Balay   if (flg) {
4140f1074feSSatish Balay     if (jac->pmax < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"P_max %G must be greater than or equal to zero",jac->pmax);
4150f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);CHKERRQ(ierr);
4160f1074feSSatish Balay   }
4170f1074feSSatish Balay 
4180f1074feSSatish Balay  ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);CHKERRQ(ierr);
4190f1074feSSatish Balay   if (flg) {
4200f1074feSSatish Balay      if (jac->agg_nl < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %G must be greater than or equal to zero",jac->agg_nl);
4210f1074feSSatish Balay 
4220f1074feSSatish Balay      ierr = HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl);CHKERRQ(ierr);
4230f1074feSSatish Balay   }
4240f1074feSSatish Balay 
4250f1074feSSatish Balay 
4260f1074feSSatish 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);
4270f1074feSSatish Balay   if (flg) {
4280f1074feSSatish Balay      if (jac->agg_num_paths < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of paths %G must be greater than or equal to 1",jac->agg_num_paths);
4290f1074feSSatish Balay 
4300f1074feSSatish Balay      ierr = HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);CHKERRQ(ierr);
4310f1074feSSatish Balay   }
4320f1074feSSatish Balay 
4330f1074feSSatish Balay 
43416d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr);
43516d9e3a6SLisandro Dalcin   if (flg) {
43616d9e3a6SLisandro Dalcin     if (jac->strongthreshold < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold);
43716d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr);
43816d9e3a6SLisandro Dalcin   }
43916d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr);
44016d9e3a6SLisandro Dalcin   if (flg) {
44116d9e3a6SLisandro Dalcin     if (jac->maxrowsum < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum);
44216d9e3a6SLisandro Dalcin     if (jac->maxrowsum > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum);
44316d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr);
44416d9e3a6SLisandro Dalcin   }
44516d9e3a6SLisandro Dalcin 
44616d9e3a6SLisandro Dalcin   /* Grid sweeps */
4470f1074feSSatish 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);
44816d9e3a6SLisandro Dalcin   if (flg) {
44916d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver,indx);CHKERRQ(ierr);
45016d9e3a6SLisandro Dalcin     /* modify the jac structure so we can view the updated options with PC_View */
45116d9e3a6SLisandro Dalcin     jac->gridsweeps[0] = indx;
4520f1074feSSatish Balay     jac->gridsweeps[1] = indx;
4530f1074feSSatish Balay     /*defaults coarse to 1 */
4540f1074feSSatish Balay     jac->gridsweeps[2] = 1;
45516d9e3a6SLisandro Dalcin   }
4560f1074feSSatish Balay 
4570f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx ,&flg);CHKERRQ(ierr);
45816d9e3a6SLisandro Dalcin   if (flg) {
45916d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 1);CHKERRQ(ierr);
4600f1074feSSatish Balay     jac->gridsweeps[0] = indx;
46116d9e3a6SLisandro Dalcin   }
46216d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr);
46316d9e3a6SLisandro Dalcin   if (flg) {
46416d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 2);CHKERRQ(ierr);
4650f1074feSSatish Balay     jac->gridsweeps[1] = indx;
46616d9e3a6SLisandro Dalcin   }
4670f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr);
46816d9e3a6SLisandro Dalcin   if (flg) {
46916d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 3);CHKERRQ(ierr);
4700f1074feSSatish Balay     jac->gridsweeps[2] = indx;
47116d9e3a6SLisandro Dalcin   }
47216d9e3a6SLisandro Dalcin 
47316d9e3a6SLisandro Dalcin   /* Relax type */
4740f1074feSSatish 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);
47516d9e3a6SLisandro Dalcin   if (flg) {
4760f1074feSSatish Balay     jac->relaxtype[0] = jac->relaxtype[1]  = indx;
47716d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, indx);CHKERRQ(ierr);
4780f1074feSSatish Balay     /* by default, coarse type set to 9 */
4790f1074feSSatish Balay     jac->relaxtype[2] = 9;
4800f1074feSSatish Balay 
48116d9e3a6SLisandro Dalcin   }
4820f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
48316d9e3a6SLisandro Dalcin   if (flg) {
48416d9e3a6SLisandro Dalcin     jac->relaxtype[0] = indx;
48516d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 1);CHKERRQ(ierr);
48616d9e3a6SLisandro Dalcin   }
4870f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
48816d9e3a6SLisandro Dalcin   if (flg) {
4890f1074feSSatish Balay     jac->relaxtype[1] = indx;
49016d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 2);CHKERRQ(ierr);
49116d9e3a6SLisandro Dalcin   }
49216d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr);
49316d9e3a6SLisandro Dalcin   if (flg) {
4940f1074feSSatish Balay     jac->relaxtype[2] = indx;
49516d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 3);CHKERRQ(ierr);
49616d9e3a6SLisandro Dalcin   }
49716d9e3a6SLisandro Dalcin 
49816d9e3a6SLisandro Dalcin   /* Relaxation Weight */
49916d9e3a6SLisandro 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);
50016d9e3a6SLisandro Dalcin   if (flg) {
50116d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetRelaxWt(jac->hsolver,tmpdbl);CHKERRQ(ierr);
50216d9e3a6SLisandro Dalcin     jac->relaxweight = tmpdbl;
50316d9e3a6SLisandro Dalcin   }
50416d9e3a6SLisandro Dalcin 
50516d9e3a6SLisandro Dalcin   n=2;
50616d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
50716d9e3a6SLisandro 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);
50816d9e3a6SLisandro Dalcin   if (flg) {
50916d9e3a6SLisandro Dalcin     if (n == 2) {
51016d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
51116d9e3a6SLisandro Dalcin       ierr = HYPRE_BoomerAMGSetLevelRelaxWt(jac->hsolver,twodbl[0],indx);CHKERRQ(ierr);
51216d9e3a6SLisandro Dalcin     } else {
51316d9e3a6SLisandro Dalcin       SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Relax weight level: you must provide 2 values separated by a comma (and no space), you provided %d",n);
51416d9e3a6SLisandro Dalcin     }
51516d9e3a6SLisandro Dalcin   }
51616d9e3a6SLisandro Dalcin 
51716d9e3a6SLisandro Dalcin   /* Outer relaxation Weight */
51816d9e3a6SLisandro 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);
51916d9e3a6SLisandro Dalcin   if (flg) {
52016d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetOuterWt( jac->hsolver, tmpdbl);CHKERRQ(ierr);
52116d9e3a6SLisandro Dalcin     jac->outerrelaxweight = tmpdbl;
52216d9e3a6SLisandro Dalcin   }
52316d9e3a6SLisandro Dalcin 
52416d9e3a6SLisandro Dalcin   n=2;
52516d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
52616d9e3a6SLisandro 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);
52716d9e3a6SLisandro Dalcin   if (flg) {
52816d9e3a6SLisandro Dalcin     if (n == 2) {
52916d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
53016d9e3a6SLisandro Dalcin       ierr = HYPRE_BoomerAMGSetLevelOuterWt( jac->hsolver, twodbl[0], indx);CHKERRQ(ierr);
53116d9e3a6SLisandro Dalcin     } else {
53216d9e3a6SLisandro Dalcin       SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Relax weight outer level: You must provide 2 values separated by a comma (and no space), you provided %d",n);
53316d9e3a6SLisandro Dalcin     }
53416d9e3a6SLisandro Dalcin   }
53516d9e3a6SLisandro Dalcin 
53616d9e3a6SLisandro Dalcin   /* the Relax Order */
53716d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
53816d9e3a6SLisandro Dalcin 
53916d9e3a6SLisandro Dalcin   if (flg) {
54016d9e3a6SLisandro Dalcin     jac->relaxorder = 0;
54116d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr);
54216d9e3a6SLisandro Dalcin   }
54316d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr);
54416d9e3a6SLisandro Dalcin   if (flg) {
54516d9e3a6SLisandro Dalcin     jac->measuretype = indx;
54616d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr);
54716d9e3a6SLisandro Dalcin   }
5480f1074feSSatish Balay   /* update list length 3/07 */
5490f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,11,HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr);
55016d9e3a6SLisandro Dalcin   if (flg) {
55116d9e3a6SLisandro Dalcin     jac->coarsentype = indx;
55216d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr);
55316d9e3a6SLisandro Dalcin   }
5540f1074feSSatish Balay 
5550f1074feSSatish Balay   /* new 3/07 */
5560f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,14,HYPREBoomerAMGInterpType[0],&indx,&flg);CHKERRQ(ierr);
5570f1074feSSatish Balay   if (flg) {
5580f1074feSSatish Balay     jac->interptype = indx;
5590f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);CHKERRQ(ierr);
5600f1074feSSatish Balay   }
5610f1074feSSatish Balay 
56290d69ab7SBarry Smith   flg  = PETSC_FALSE;
56390d69ab7SBarry Smith   ierr = PetscOptionsTruth("-pc_hypre_boomeramg_print_statistics","Print statistics","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
56416d9e3a6SLisandro Dalcin   if (flg) {
56516d9e3a6SLisandro Dalcin     int level=3;
56616d9e3a6SLisandro Dalcin     jac->printstatistics = PETSC_TRUE;
56716d9e3a6SLisandro Dalcin     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
56816d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetPrintLevel(jac->hsolver,level);CHKERRQ(ierr);
5692ae77aedSBarry Smith   }
5702ae77aedSBarry Smith 
57190d69ab7SBarry Smith   flg  = PETSC_FALSE;
57290d69ab7SBarry Smith   ierr = PetscOptionsTruth("-pc_hypre_boomeramg_print_debug","Print debug information","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
5732ae77aedSBarry Smith   if (flg) {
5742ae77aedSBarry Smith     int level=3;
5752ae77aedSBarry Smith     jac->printstatistics = PETSC_TRUE;
5762ae77aedSBarry Smith     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
57716d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetDebugFlag(jac->hsolver,level);CHKERRQ(ierr);
57816d9e3a6SLisandro Dalcin   }
5798f87f92bSBarry Smith 
5808f87f92bSBarry Smith   ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5818f87f92bSBarry Smith   if (flg && tmp_truth) {
5828f87f92bSBarry Smith     jac->nodal_coarsen = 1;
5838f87f92bSBarry Smith     ierr = HYPRE_BoomerAMGSetNodal(jac->hsolver,1);CHKERRQ(ierr);
5848f87f92bSBarry Smith   }
5858f87f92bSBarry Smith 
5868f87f92bSBarry Smith   ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5878f87f92bSBarry Smith   if (flg && tmp_truth) {
5888f87f92bSBarry Smith     PetscInt tmp_int;
5898f87f92bSBarry Smith     ierr = PetscOptionsInt( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);CHKERRQ(ierr);
5908f87f92bSBarry Smith     if (flg) jac->nodal_relax_levels = tmp_int;
5918f87f92bSBarry Smith     ierr = HYPRE_BoomerAMGSetSmoothType(jac->hsolver,6);CHKERRQ(ierr);
5928f87f92bSBarry Smith     ierr = HYPRE_BoomerAMGSetDomainType(jac->hsolver,1);CHKERRQ(ierr);
5938f87f92bSBarry Smith     ierr = HYPRE_BoomerAMGSetOverlap(jac->hsolver,0);CHKERRQ(ierr);
5948f87f92bSBarry Smith     ierr = HYPRE_BoomerAMGSetSmoothNumLevels(jac->hsolver,jac->nodal_relax_levels);CHKERRQ(ierr);
5958f87f92bSBarry Smith   }
5968f87f92bSBarry Smith 
59716d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
59816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
59916d9e3a6SLisandro Dalcin }
60016d9e3a6SLisandro Dalcin 
60116d9e3a6SLisandro Dalcin #undef __FUNCT__
60216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG"
6034d0a8057SBarry Smith static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscInt *outits,PCRichardsonConvergedReason *reason)
60416d9e3a6SLisandro Dalcin {
60516d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
60616d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
6074d0a8057SBarry Smith   int            oits;
60816d9e3a6SLisandro Dalcin 
60916d9e3a6SLisandro Dalcin   PetscFunctionBegin;
61016d9e3a6SLisandro Dalcin   ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,its*jac->maxiter);CHKERRQ(ierr);
6114d0a8057SBarry Smith   ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,rtol);CHKERRQ(ierr);
61216d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_TRUE;
61316d9e3a6SLisandro Dalcin   ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr);
61416d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_FALSE;
6154d0a8057SBarry Smith   ierr = HYPRE_BoomerAMGGetNumIterations(jac->hsolver,&oits);CHKERRQ(ierr);
6164d0a8057SBarry Smith   *outits = oits;
6174d0a8057SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
6184d0a8057SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
61916d9e3a6SLisandro Dalcin   ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
62016d9e3a6SLisandro Dalcin   ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
62116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
62216d9e3a6SLisandro Dalcin }
62316d9e3a6SLisandro Dalcin 
62416d9e3a6SLisandro Dalcin 
62516d9e3a6SLisandro Dalcin #undef __FUNCT__
62616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_BoomerAMG"
62716d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
62816d9e3a6SLisandro Dalcin {
62916d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
63016d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
63116d9e3a6SLisandro Dalcin   PetscTruth     iascii;
63216d9e3a6SLisandro Dalcin 
63316d9e3a6SLisandro Dalcin   PetscFunctionBegin;
63416d9e3a6SLisandro Dalcin   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
63516d9e3a6SLisandro Dalcin   if (iascii) {
63616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr);
63716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr);
63816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr);
63916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr);
64016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr);
64116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr);
6420f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation truncation factor %G\n",jac->truncfactor);CHKERRQ(ierr);
6430f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);CHKERRQ(ierr);
6440f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);CHKERRQ(ierr);
6450f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);CHKERRQ(ierr);
6460f1074feSSatish Balay 
64716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr);
64816d9e3a6SLisandro Dalcin 
6490f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps down         %d\n",jac->gridsweeps[0]);CHKERRQ(ierr);
6500f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps up           %d\n",jac->gridsweeps[1]);CHKERRQ(ierr);
6510f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on coarse    %d\n",jac->gridsweeps[2]);CHKERRQ(ierr);
65216d9e3a6SLisandro Dalcin 
6530f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr);
6540f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr);
6550f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr);
65616d9e3a6SLisandro Dalcin 
65716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax weight  (all)      %G\n",jac->relaxweight);CHKERRQ(ierr);
65816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr);
65916d9e3a6SLisandro Dalcin 
66016d9e3a6SLisandro Dalcin     if (jac->relaxorder) {
66116d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr);
66216d9e3a6SLisandro Dalcin     } else {
66316d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr);
66416d9e3a6SLisandro Dalcin     }
66516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr);
66616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr);
6670f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation type  %s\n",HYPREBoomerAMGInterpType[jac->interptype]);CHKERRQ(ierr);
6688f87f92bSBarry Smith     if (jac->nodal_coarsen) {
6698f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal())\n");CHKERRQ(ierr);
6708f87f92bSBarry Smith     }
6718f87f92bSBarry Smith     if (jac->nodal_relax) {
6728f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);CHKERRQ(ierr);
6738f87f92bSBarry Smith     }
67416d9e3a6SLisandro Dalcin   }
67516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
67616d9e3a6SLisandro Dalcin }
67716d9e3a6SLisandro Dalcin 
67816d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
67916d9e3a6SLisandro Dalcin #undef __FUNCT__
68016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails"
68116d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc)
68216d9e3a6SLisandro Dalcin {
68316d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
68416d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
68516d9e3a6SLisandro Dalcin   int            indx;
68616d9e3a6SLisandro Dalcin   PetscTruth     flag;
68716d9e3a6SLisandro Dalcin   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};
68816d9e3a6SLisandro Dalcin 
68916d9e3a6SLisandro Dalcin   PetscFunctionBegin;
69016d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr);
69116d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr);
69216d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr);
69316d9e3a6SLisandro Dalcin   if (flag) {
69416d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr);
69516d9e3a6SLisandro Dalcin   }
69616d9e3a6SLisandro Dalcin 
69716d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr);
69816d9e3a6SLisandro Dalcin   if (flag) {
69916d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr);
70016d9e3a6SLisandro Dalcin   }
70116d9e3a6SLisandro Dalcin 
70216d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr);
70316d9e3a6SLisandro Dalcin   if (flag) {
70416d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr);
70516d9e3a6SLisandro Dalcin   }
70616d9e3a6SLisandro Dalcin 
70716d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,&flag);CHKERRQ(ierr);
70816d9e3a6SLisandro Dalcin   if (flag) {
70916d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr);
71016d9e3a6SLisandro Dalcin   }
71116d9e3a6SLisandro Dalcin 
71216d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,&flag);CHKERRQ(ierr);
71316d9e3a6SLisandro Dalcin   if (flag) {
71416d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr);
71516d9e3a6SLisandro Dalcin   }
71616d9e3a6SLisandro Dalcin 
71716d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);CHKERRQ(ierr);
71816d9e3a6SLisandro Dalcin   if (flag) {
71916d9e3a6SLisandro Dalcin     jac->symt = indx;
72016d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr);
72116d9e3a6SLisandro Dalcin   }
72216d9e3a6SLisandro Dalcin 
72316d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
72416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
72516d9e3a6SLisandro Dalcin }
72616d9e3a6SLisandro Dalcin 
72716d9e3a6SLisandro Dalcin #undef __FUNCT__
72816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_ParaSails"
72916d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
73016d9e3a6SLisandro Dalcin {
73116d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
73216d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
73316d9e3a6SLisandro Dalcin   PetscTruth     iascii;
73416d9e3a6SLisandro Dalcin   const char     *symt = 0;;
73516d9e3a6SLisandro Dalcin 
73616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
73716d9e3a6SLisandro Dalcin   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
73816d9e3a6SLisandro Dalcin   if (iascii) {
73916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");CHKERRQ(ierr);
74016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr);
74116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr);
74216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr);
74316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr);
74416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscTruths[jac->ruse]);CHKERRQ(ierr);
74516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscTruths[jac->logging]);CHKERRQ(ierr);
74616d9e3a6SLisandro Dalcin     if (!jac->symt) {
74716d9e3a6SLisandro Dalcin       symt = "nonsymmetric matrix and preconditioner";
74816d9e3a6SLisandro Dalcin     } else if (jac->symt == 1) {
74916d9e3a6SLisandro Dalcin       symt = "SPD matrix and preconditioner";
75016d9e3a6SLisandro Dalcin     } else if (jac->symt == 2) {
75116d9e3a6SLisandro Dalcin       symt = "nonsymmetric matrix but SPD preconditioner";
75216d9e3a6SLisandro Dalcin     } else {
75316d9e3a6SLisandro Dalcin       SETERRQ1(PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
75416d9e3a6SLisandro Dalcin     }
75516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr);
75616d9e3a6SLisandro Dalcin   }
75716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
75816d9e3a6SLisandro Dalcin }
75916d9e3a6SLisandro Dalcin /* ---------------------------------------------------------------------------------*/
76016d9e3a6SLisandro Dalcin 
76116d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
76216d9e3a6SLisandro Dalcin #undef __FUNCT__
76316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType_HYPRE"
76416d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType_HYPRE(PC pc,const char *name[])
76516d9e3a6SLisandro Dalcin {
76616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
76716d9e3a6SLisandro Dalcin 
76816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
76916d9e3a6SLisandro Dalcin   *name = jac->hypre_type;
77016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
77116d9e3a6SLisandro Dalcin }
77216d9e3a6SLisandro Dalcin EXTERN_C_END
77316d9e3a6SLisandro Dalcin 
77416d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
77516d9e3a6SLisandro Dalcin #undef __FUNCT__
77616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType_HYPRE"
77716d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType_HYPRE(PC pc,const char name[])
77816d9e3a6SLisandro Dalcin {
77916d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
78016d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
78116d9e3a6SLisandro Dalcin   PetscTruth     flag;
78216d9e3a6SLisandro Dalcin 
78316d9e3a6SLisandro Dalcin   PetscFunctionBegin;
78416d9e3a6SLisandro Dalcin   if (jac->hypre_type) {
78516d9e3a6SLisandro Dalcin     ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr);
78616d9e3a6SLisandro Dalcin     if (!flag) SETERRQ(PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
78716d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
78816d9e3a6SLisandro Dalcin   } else {
78916d9e3a6SLisandro Dalcin     ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr);
79016d9e3a6SLisandro Dalcin   }
79116d9e3a6SLisandro Dalcin 
79216d9e3a6SLisandro Dalcin   jac->maxiter            = PETSC_DEFAULT;
79316d9e3a6SLisandro Dalcin   jac->tol                = PETSC_DEFAULT;
79416d9e3a6SLisandro Dalcin   jac->printstatistics    = PetscLogPrintInfo;
79516d9e3a6SLisandro Dalcin 
79616d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr);
79716d9e3a6SLisandro Dalcin   if (flag) {
79816d9e3a6SLisandro Dalcin     ierr                    = HYPRE_ParCSRPilutCreate(jac->comm_hypre,&jac->hsolver);
79916d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
80016d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Pilut;
80116d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParCSRPilutDestroy;
80216d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParCSRPilutSetup;
80316d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParCSRPilutSolve;
80416d9e3a6SLisandro Dalcin     jac->factorrowsize      = PETSC_DEFAULT;
80516d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
80616d9e3a6SLisandro Dalcin   }
80716d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr);
80816d9e3a6SLisandro Dalcin   if (flag) {
80916d9e3a6SLisandro Dalcin     ierr                    = HYPRE_ParaSailsCreate(jac->comm_hypre,&jac->hsolver);
81016d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
81116d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_ParaSails;
81216d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParaSailsDestroy;
81316d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParaSailsSetup;
81416d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParaSailsSolve;
81516d9e3a6SLisandro Dalcin     /* initialize */
81616d9e3a6SLisandro Dalcin     jac->nlevels     = 1;
81716d9e3a6SLisandro Dalcin     jac->threshhold  = .1;
81816d9e3a6SLisandro Dalcin     jac->filter      = .1;
81916d9e3a6SLisandro Dalcin     jac->loadbal     = 0;
82016d9e3a6SLisandro Dalcin     if (PetscLogPrintInfo) {
82116d9e3a6SLisandro Dalcin       jac->logging     = (int) PETSC_TRUE;
82216d9e3a6SLisandro Dalcin     } else {
82316d9e3a6SLisandro Dalcin       jac->logging     = (int) PETSC_FALSE;
82416d9e3a6SLisandro Dalcin     }
82516d9e3a6SLisandro Dalcin     jac->ruse = (int) PETSC_FALSE;
82616d9e3a6SLisandro Dalcin     jac->symt   = 0;
82716d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr);
82816d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr);
82916d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr);
83016d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr);
83116d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr);
83216d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr);
83316d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
83416d9e3a6SLisandro Dalcin   }
83516d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr);
83616d9e3a6SLisandro Dalcin   if (flag) {
83716d9e3a6SLisandro Dalcin     ierr                    = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver);
83816d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
83916d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Euclid;
84016d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_EuclidDestroy;
84116d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_EuclidSetup;
84216d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_EuclidSolve;
84316d9e3a6SLisandro Dalcin     /* initialization */
84416d9e3a6SLisandro Dalcin     jac->bjilu              = PETSC_FALSE;
84516d9e3a6SLisandro Dalcin     jac->levels             = 1;
84616d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
84716d9e3a6SLisandro Dalcin   }
84816d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr);
84916d9e3a6SLisandro Dalcin   if (flag) {
85016d9e3a6SLisandro Dalcin     ierr                     = HYPRE_BoomerAMGCreate(&jac->hsolver);
85116d9e3a6SLisandro Dalcin     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
85216d9e3a6SLisandro Dalcin     pc->ops->view            = PCView_HYPRE_BoomerAMG;
85316d9e3a6SLisandro Dalcin     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
85416d9e3a6SLisandro Dalcin     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
85516d9e3a6SLisandro Dalcin     jac->destroy             = HYPRE_BoomerAMGDestroy;
85616d9e3a6SLisandro Dalcin     jac->setup               = HYPRE_BoomerAMGSetup;
85716d9e3a6SLisandro Dalcin     jac->solve               = HYPRE_BoomerAMGSolve;
85816d9e3a6SLisandro Dalcin     jac->applyrichardson     = PETSC_FALSE;
85916d9e3a6SLisandro Dalcin     /* these defaults match the hypre defaults */
86016d9e3a6SLisandro Dalcin     jac->cycletype        = 1;
86116d9e3a6SLisandro Dalcin     jac->maxlevels        = 25;
86216d9e3a6SLisandro Dalcin     jac->maxiter          = 1;
8638f87f92bSBarry Smith     jac->tol              = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */
86416d9e3a6SLisandro Dalcin     jac->truncfactor      = 0.0;
86516d9e3a6SLisandro Dalcin     jac->strongthreshold  = .25;
86616d9e3a6SLisandro Dalcin     jac->maxrowsum        = .9;
86716d9e3a6SLisandro Dalcin     jac->coarsentype      = 6;
86816d9e3a6SLisandro Dalcin     jac->measuretype      = 0;
8690f1074feSSatish Balay     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = 1;
8708f87f92bSBarry Smith     jac->relaxtype[0]     = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */
8710f1074feSSatish Balay     jac->relaxtype[2]     = 9; /*G.E. */
87216d9e3a6SLisandro Dalcin     jac->relaxweight      = 1.0;
87316d9e3a6SLisandro Dalcin     jac->outerrelaxweight = 1.0;
87416d9e3a6SLisandro Dalcin     jac->relaxorder       = 1;
8750f1074feSSatish Balay     jac->interptype       = 0;
8760f1074feSSatish Balay     jac->agg_nl           = 0;
8770f1074feSSatish Balay     jac->pmax             = 0;
8780f1074feSSatish Balay     jac->truncfactor      = 0.0;
8790f1074feSSatish Balay     jac->agg_num_paths    = 1;
8808f87f92bSBarry Smith 
8818f87f92bSBarry Smith     jac->nodal_coarsen    = 0;
8828f87f92bSBarry Smith     jac->nodal_relax      = PETSC_FALSE;
8838f87f92bSBarry Smith     jac->nodal_relax_levels = 1;
88416d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr);
88516d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr);
88616d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
88716d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
88816d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr);
88916d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr);
89016d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr);
89116d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr);
89216d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr);
89316d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr);
8940f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);CHKERRQ(ierr);
8950f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl);
8960f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);CHKERRQ(ierr);
8970f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);CHKERRQ(ierr);
8980f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, jac->relaxtype[0]);  /*defaults coarse to 9*/
8990f1074feSSatish Balay     ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver, jac->gridsweeps[0]); /*defaults coarse to 1 */
90016d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
90116d9e3a6SLisandro Dalcin   }
90216d9e3a6SLisandro Dalcin   ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr);
90316d9e3a6SLisandro Dalcin   jac->hypre_type = PETSC_NULL;
90416d9e3a6SLisandro Dalcin   SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
90516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
90616d9e3a6SLisandro Dalcin }
90716d9e3a6SLisandro Dalcin EXTERN_C_END
90816d9e3a6SLisandro Dalcin 
90916d9e3a6SLisandro Dalcin /*
91016d9e3a6SLisandro Dalcin     It only gets here if the HYPRE type has not been set before the call to
91116d9e3a6SLisandro Dalcin    ...SetFromOptions() which actually is most of the time
91216d9e3a6SLisandro Dalcin */
91316d9e3a6SLisandro Dalcin #undef __FUNCT__
91416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE"
91516d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE(PC pc)
91616d9e3a6SLisandro Dalcin {
91716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
91816d9e3a6SLisandro Dalcin   int            indx;
91916d9e3a6SLisandro Dalcin   const char     *type[] = {"pilut","parasails","boomeramg","euclid"};
92016d9e3a6SLisandro Dalcin   PetscTruth     flg;
92116d9e3a6SLisandro Dalcin 
92216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
92316d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr);
92402a17cd4SBarry Smith   ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"boomeramg",&indx,&flg);CHKERRQ(ierr);
92516d9e3a6SLisandro Dalcin   if (flg) {
92616d9e3a6SLisandro Dalcin     ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr);
92702a17cd4SBarry Smith   } else {
92802a17cd4SBarry Smith     ierr = PCHYPRESetType_HYPRE(pc,"boomeramg");CHKERRQ(ierr);
92916d9e3a6SLisandro Dalcin   }
93016d9e3a6SLisandro Dalcin   if (pc->ops->setfromoptions) {
93116d9e3a6SLisandro Dalcin     ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr);
93216d9e3a6SLisandro Dalcin   }
93316d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
93416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
93516d9e3a6SLisandro Dalcin }
93616d9e3a6SLisandro Dalcin 
93716d9e3a6SLisandro Dalcin #undef __FUNCT__
93816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType"
93916d9e3a6SLisandro Dalcin /*@C
94016d9e3a6SLisandro Dalcin      PCHYPRESetType - Sets which hypre preconditioner you wish to use
94116d9e3a6SLisandro Dalcin 
94216d9e3a6SLisandro Dalcin    Input Parameters:
94316d9e3a6SLisandro Dalcin +     pc - the preconditioner context
94416d9e3a6SLisandro Dalcin -     name - either  pilut, parasails, boomeramg, euclid
94516d9e3a6SLisandro Dalcin 
94616d9e3a6SLisandro Dalcin    Options Database Keys:
94716d9e3a6SLisandro Dalcin    -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
94816d9e3a6SLisandro Dalcin 
94916d9e3a6SLisandro Dalcin    Level: intermediate
95016d9e3a6SLisandro Dalcin 
95116d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
95216d9e3a6SLisandro Dalcin            PCHYPRE
95316d9e3a6SLisandro Dalcin 
95416d9e3a6SLisandro Dalcin @*/
95516d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC pc,const char name[])
95616d9e3a6SLisandro Dalcin {
95716d9e3a6SLisandro Dalcin   PetscErrorCode ierr,(*f)(PC,const char[]);
95816d9e3a6SLisandro Dalcin 
95916d9e3a6SLisandro Dalcin   PetscFunctionBegin;
96016d9e3a6SLisandro Dalcin   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
96116d9e3a6SLisandro Dalcin   PetscValidCharPointer(name,2);
96216d9e3a6SLisandro Dalcin   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPRESetType_C",(void (**)(void))&f);CHKERRQ(ierr);
96316d9e3a6SLisandro Dalcin   if (f) {
96416d9e3a6SLisandro Dalcin     ierr = (*f)(pc,name);CHKERRQ(ierr);
96516d9e3a6SLisandro Dalcin   }
96616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
96716d9e3a6SLisandro Dalcin }
96816d9e3a6SLisandro Dalcin 
96916d9e3a6SLisandro Dalcin #undef __FUNCT__
97016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType"
97116d9e3a6SLisandro Dalcin /*@C
97216d9e3a6SLisandro Dalcin      PCHYPREGetType - Gets which hypre preconditioner you are using
97316d9e3a6SLisandro Dalcin 
97416d9e3a6SLisandro Dalcin    Input Parameter:
97516d9e3a6SLisandro Dalcin .     pc - the preconditioner context
97616d9e3a6SLisandro Dalcin 
97716d9e3a6SLisandro Dalcin    Output Parameter:
97816d9e3a6SLisandro Dalcin .     name - either  pilut, parasails, boomeramg, euclid
97916d9e3a6SLisandro Dalcin 
98016d9e3a6SLisandro Dalcin    Level: intermediate
98116d9e3a6SLisandro Dalcin 
98216d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC,
98316d9e3a6SLisandro Dalcin            PCHYPRE
98416d9e3a6SLisandro Dalcin 
98516d9e3a6SLisandro Dalcin @*/
98616d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC pc,const char *name[])
98716d9e3a6SLisandro Dalcin {
98816d9e3a6SLisandro Dalcin   PetscErrorCode ierr,(*f)(PC,const char*[]);
98916d9e3a6SLisandro Dalcin 
99016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
99116d9e3a6SLisandro Dalcin   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
99216d9e3a6SLisandro Dalcin   PetscValidPointer(name,2);
99316d9e3a6SLisandro Dalcin   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPREGetType_C",(void (**)(void))&f);CHKERRQ(ierr);
99416d9e3a6SLisandro Dalcin   if (f) {
99516d9e3a6SLisandro Dalcin     ierr = (*f)(pc,name);CHKERRQ(ierr);
99616d9e3a6SLisandro Dalcin   }
99716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
99816d9e3a6SLisandro Dalcin }
99916d9e3a6SLisandro Dalcin 
100016d9e3a6SLisandro Dalcin /*MC
100116d9e3a6SLisandro Dalcin      PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre
100216d9e3a6SLisandro Dalcin 
100316d9e3a6SLisandro Dalcin    Options Database Keys:
100416d9e3a6SLisandro Dalcin +   -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
100516d9e3a6SLisandro Dalcin -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
100616d9e3a6SLisandro Dalcin           preconditioner
100716d9e3a6SLisandro Dalcin 
100816d9e3a6SLisandro Dalcin    Level: intermediate
100916d9e3a6SLisandro Dalcin 
101016d9e3a6SLisandro Dalcin    Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()),
101116d9e3a6SLisandro Dalcin           the many hypre options can ONLY be set via the options database (e.g. the command line
101216d9e3a6SLisandro Dalcin           or with PetscOptionsSetValue(), there are no functions to set them)
101316d9e3a6SLisandro Dalcin 
101416d9e3a6SLisandro Dalcin           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations
10150f1074feSSatish Balay           (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if
10160f1074feSSatish Balay           -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner
10170f1074feSSatish Balay           (-pc_hypre_boomeramg_rtol should be set to 0.0 - the default - to strictly use a fixed number of
10188f87f92bSBarry Smith           iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations
10190f1074feSSatish Balay           and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10
10200f1074feSSatish Balay           then AT MOST twenty V-cycles of boomeramg will be called.
102116d9e3a6SLisandro Dalcin 
10220f1074feSSatish Balay            Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation
10230f1074feSSatish Balay            (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry.
10240f1074feSSatish Balay            Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi.
102516d9e3a6SLisandro Dalcin           If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly
102616d9e3a6SLisandro Dalcin           and use -ksp_max_it to control the number of V-cycles.
102716d9e3a6SLisandro Dalcin           (see the PETSc FAQ.html at the PETSc website under the Documentation tab).
102816d9e3a6SLisandro Dalcin 
102916d9e3a6SLisandro Dalcin 	  2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option
103016d9e3a6SLisandro Dalcin 	  -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L.
103116d9e3a6SLisandro Dalcin 
10329e5bc791SBarry Smith           See PCPFMG for access to the hypre Struct PFMG solver
10339e5bc791SBarry Smith 
103416d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
10359e5bc791SBarry Smith            PCHYPRESetType(), PCPFMG
103616d9e3a6SLisandro Dalcin 
103716d9e3a6SLisandro Dalcin M*/
103816d9e3a6SLisandro Dalcin 
103916d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
104016d9e3a6SLisandro Dalcin #undef __FUNCT__
104116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCCreate_HYPRE"
104216d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_HYPRE(PC pc)
104316d9e3a6SLisandro Dalcin {
104416d9e3a6SLisandro Dalcin   PC_HYPRE       *jac;
104516d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
104616d9e3a6SLisandro Dalcin 
104716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
104838f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_HYPRE,&jac);CHKERRQ(ierr);
104916d9e3a6SLisandro Dalcin   pc->data                 = jac;
105016d9e3a6SLisandro Dalcin   pc->ops->destroy         = PCDestroy_HYPRE;
105116d9e3a6SLisandro Dalcin   pc->ops->setfromoptions  = PCSetFromOptions_HYPRE;
105216d9e3a6SLisandro Dalcin   pc->ops->setup           = PCSetUp_HYPRE;
105316d9e3a6SLisandro Dalcin   pc->ops->apply           = PCApply_HYPRE;
105416d9e3a6SLisandro Dalcin   jac->comm_hypre          = MPI_COMM_NULL;
105516d9e3a6SLisandro Dalcin   jac->hypre_type          = PETSC_NULL;
105616d9e3a6SLisandro Dalcin   /* duplicate communicator for hypre */
10577adad957SLisandro Dalcin   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(jac->comm_hypre));CHKERRQ(ierr);
105816d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr);
105916d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr);
106016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
106116d9e3a6SLisandro Dalcin }
106216d9e3a6SLisandro Dalcin EXTERN_C_END
1063ebc551c0SBarry Smith 
1064f91d8e95SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------*/
1065f91d8e95SBarry Smith 
106668326731SBarry Smith /* we know we are working with a HYPRE_StructMatrix */
106768326731SBarry Smith #include "../src/mat/impls/hypre/mhyp.h"
106868326731SBarry Smith #include "private/matimpl.h"
106968326731SBarry Smith #include "private/pcimpl.h"
1070ebc551c0SBarry Smith 
1071ebc551c0SBarry Smith typedef struct {
107268326731SBarry Smith   MPI_Comm            hcomm;       /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */
1073f91d8e95SBarry Smith   HYPRE_StructSolver  hsolver;
10749e5bc791SBarry Smith 
10759e5bc791SBarry Smith   /* keep copy of PFMG options used so may view them */
10769e5bc791SBarry Smith   int                 its;
10779e5bc791SBarry Smith   double              tol;
10789e5bc791SBarry Smith   int                 relax_type;
10799e5bc791SBarry Smith   int                 rap_type;
10809e5bc791SBarry Smith   int                 num_pre_relax,num_post_relax;
1081ebc551c0SBarry Smith } PC_PFMG;
1082ebc551c0SBarry Smith 
1083ebc551c0SBarry Smith #undef __FUNCT__
1084ebc551c0SBarry Smith #define __FUNCT__ "PCSetUp_PFMG"
1085ebc551c0SBarry Smith PetscErrorCode PCSetUp_PFMG(PC pc)
1086ebc551c0SBarry Smith {
1087ebc551c0SBarry Smith   PetscErrorCode  ierr;
1088ebc551c0SBarry Smith   PC_PFMG         *ex = (PC_PFMG*) pc->data;
108968326731SBarry Smith   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
109068326731SBarry Smith   PetscTruth      flg;
1091ebc551c0SBarry Smith 
1092ebc551c0SBarry Smith   PetscFunctionBegin;
109368326731SBarry Smith   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATHYPRESTRUCT,&flg);CHKERRQ(ierr);
109468326731SBarry Smith   if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESTRUCT with this preconditioner");
1095f91d8e95SBarry Smith 
1096f91d8e95SBarry Smith   /* create the hypre solver object and set its information */
1097*c7ef6ce2SGlenn Hammond   if (ex->hsolver) {
1098*c7ef6ce2SGlenn Hammond     ierr = HYPRE_StructPFMGDestroy(ex->hsolver);CHKERRQ(ierr);
1099*c7ef6ce2SGlenn Hammond   }
1100*c7ef6ce2SGlenn Hammond   ierr = HYPRE_StructPFMGCreate(ex->hcomm,&ex->hsolver);CHKERRQ(ierr);
110168326731SBarry Smith   ierr = HYPRE_StructPFMGSetup(ex->hsolver,mx->hmat,mx->hb,mx->hx);CHKERRQ(ierr);
1102a0324ebeSBarry Smith   ierr = HYPRE_StructPFMGSetZeroGuess(ex->hsolver);CHKERRQ(ierr);
1103a0324ebeSBarry Smith 
1104ebc551c0SBarry Smith   PetscFunctionReturn(0);
1105ebc551c0SBarry Smith }
1106ebc551c0SBarry Smith 
1107ebc551c0SBarry Smith #undef __FUNCT__
1108ebc551c0SBarry Smith #define __FUNCT__ "PCDestroy_PFMG"
1109ebc551c0SBarry Smith PetscErrorCode PCDestroy_PFMG(PC pc)
1110ebc551c0SBarry Smith {
1111ebc551c0SBarry Smith   PetscErrorCode ierr;
1112f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1113ebc551c0SBarry Smith 
1114ebc551c0SBarry Smith   PetscFunctionBegin;
1115f91d8e95SBarry Smith   if (ex->hsolver) {ierr = HYPRE_StructPFMGDestroy(ex->hsolver);CHKERRQ(ierr);}
1116f91d8e95SBarry Smith   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1117f91d8e95SBarry Smith   ierr = PetscFree(ex);CHKERRQ(ierr);
1118ebc551c0SBarry Smith   PetscFunctionReturn(0);
1119ebc551c0SBarry Smith }
1120ebc551c0SBarry Smith 
11219e5bc791SBarry Smith static const char *PFMGRelaxType[]   = {"Jacobi","Weighted-Jacobi","symmetric-Red/Black-Gauss-Seidel","Red/Black-Gauss-Seidel"};
11229e5bc791SBarry Smith static const char *PFMGRAPType[]   = {"Galerkin","non-Galerkin"};
11239e5bc791SBarry Smith 
1124ebc551c0SBarry Smith #undef __FUNCT__
1125ebc551c0SBarry Smith #define __FUNCT__ "PCView_PFMG"
1126ebc551c0SBarry Smith PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer)
1127ebc551c0SBarry Smith {
1128ebc551c0SBarry Smith   PetscErrorCode ierr;
1129ebc551c0SBarry Smith   PetscTruth     iascii;
1130f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1131ebc551c0SBarry Smith 
1132ebc551c0SBarry Smith   PetscFunctionBegin;
1133ebc551c0SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
11349e5bc791SBarry Smith   if (iascii) {
11359e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG preconditioning\n");CHKERRQ(ierr);
11369e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
11379e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
11389e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
11399e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: RAP type %s\n",PFMGRAPType[ex->rap_type]);CHKERRQ(ierr);
11409e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
11419e5bc791SBarry Smith   }
1142ebc551c0SBarry Smith   PetscFunctionReturn(0);
1143ebc551c0SBarry Smith }
1144ebc551c0SBarry Smith 
11459e5bc791SBarry Smith 
1146ebc551c0SBarry Smith #undef __FUNCT__
1147ebc551c0SBarry Smith #define __FUNCT__ "PCSetFromOptions_PFMG"
1148ebc551c0SBarry Smith PetscErrorCode PCSetFromOptions_PFMG(PC pc)
1149ebc551c0SBarry Smith {
1150ebc551c0SBarry Smith   PetscErrorCode ierr;
1151f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
11529e5bc791SBarry Smith   PetscTruth     flg;
1153ebc551c0SBarry Smith 
1154ebc551c0SBarry Smith   PetscFunctionBegin;
1155ebc551c0SBarry Smith   ierr = PetscOptionsHead("PFMG options");CHKERRQ(ierr);
11569e5bc791SBarry Smith   ierr = PetscOptionsTruth("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
115768326731SBarry Smith   if (flg) {
1158a0324ebeSBarry Smith     int level=3;
115968326731SBarry Smith     ierr = HYPRE_StructPFMGSetPrintLevel(ex->hsolver,level);CHKERRQ(ierr);
116068326731SBarry Smith   }
11619e5bc791SBarry 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);
11629e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetMaxIter(ex->hsolver,ex->its);CHKERRQ(ierr);
11639e5bc791SBarry 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);
11649e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetNumPreRelax(ex->hsolver,ex->num_pre_relax);CHKERRQ(ierr);
11659e5bc791SBarry 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);
11669e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetNumPostRelax(ex->hsolver,ex->num_post_relax);CHKERRQ(ierr);
11679e5bc791SBarry Smith 
11689e5bc791SBarry Smith   ierr = PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr);
11699e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetTol(ex->hsolver,ex->tol);CHKERRQ(ierr);
11709e5bc791SBarry 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);
11719e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetRelaxType(ex->hsolver, ex->relax_type);CHKERRQ(ierr);
11729e5bc791SBarry Smith   ierr = PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,2,PFMGRAPType[ex->rap_type],&ex->rap_type,PETSC_NULL);CHKERRQ(ierr);
11739e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetRAPType(ex->hsolver, ex->rap_type);CHKERRQ(ierr);
1174ebc551c0SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1175ebc551c0SBarry Smith   PetscFunctionReturn(0);
1176ebc551c0SBarry Smith }
1177ebc551c0SBarry Smith 
1178f91d8e95SBarry Smith #undef __FUNCT__
1179f91d8e95SBarry Smith #define __FUNCT__ "PCApply_PFMG"
1180f91d8e95SBarry Smith PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y)
1181f91d8e95SBarry Smith {
1182f91d8e95SBarry Smith   PetscErrorCode  ierr;
1183f91d8e95SBarry Smith   PC_PFMG         *ex = (PC_PFMG*) pc->data;
1184f91d8e95SBarry Smith   PetscScalar     *xx,*yy;
1185f91d8e95SBarry Smith   int             ilower[3],iupper[3];
118668326731SBarry Smith   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
1187f91d8e95SBarry Smith 
1188f91d8e95SBarry Smith   PetscFunctionBegin;
118968326731SBarry Smith   ierr = DAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1190f91d8e95SBarry Smith   iupper[0] += ilower[0] - 1;
1191f91d8e95SBarry Smith   iupper[1] += ilower[1] - 1;
1192f91d8e95SBarry Smith   iupper[2] += ilower[2] - 1;
1193f91d8e95SBarry Smith 
1194f91d8e95SBarry Smith   /* copy x values over to hypre */
1195*c7ef6ce2SGlenn Hammond   ierr = HYPRE_StructVectorSetConstantValues(mx->hb,0.0);CHKERRQ(ierr);
1196f91d8e95SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
119768326731SBarry Smith   ierr = HYPRE_StructVectorSetBoxValues(mx->hb,ilower,iupper,xx);CHKERRQ(ierr);
1198f91d8e95SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
119968326731SBarry Smith   ierr = HYPRE_StructVectorAssemble(mx->hb);CHKERRQ(ierr);
1200f91d8e95SBarry Smith 
120168326731SBarry Smith   ierr = HYPRE_StructPFMGSolve(ex->hsolver,mx->hmat,mx->hb,mx->hx);CHKERRQ(ierr);
1202f91d8e95SBarry Smith 
1203f91d8e95SBarry Smith   /* copy solution values back to PETSc */
1204f91d8e95SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1205a0324ebeSBarry Smith   ierr = HYPRE_StructVectorGetBoxValues(mx->hx,ilower,iupper,yy);CHKERRQ(ierr);
1206f91d8e95SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1207f91d8e95SBarry Smith   PetscFunctionReturn(0);
1208f91d8e95SBarry Smith }
1209f91d8e95SBarry Smith 
12109e5bc791SBarry Smith #undef __FUNCT__
12119e5bc791SBarry Smith #define __FUNCT__ "PCApplyRichardson_PFMG"
12129e5bc791SBarry Smith static PetscErrorCode PCApplyRichardson_PFMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscInt *outits,PCRichardsonConvergedReason *reason)
12139e5bc791SBarry Smith {
12149e5bc791SBarry Smith   PC_PFMG        *jac = (PC_PFMG*)pc->data;
12159e5bc791SBarry Smith   PetscErrorCode ierr;
12169e5bc791SBarry Smith   int            oits;
12179e5bc791SBarry Smith 
12189e5bc791SBarry Smith   PetscFunctionBegin;
12199e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetMaxIter(jac->hsolver,its*jac->its);CHKERRQ(ierr);
12209e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetTol(jac->hsolver,rtol);CHKERRQ(ierr);
12219e5bc791SBarry Smith 
12229e5bc791SBarry Smith   ierr = PCApply_PFMG(pc,b,y);CHKERRQ(ierr);
12239e5bc791SBarry Smith   ierr = HYPRE_StructPFMGGetNumIterations(jac->hsolver,&oits);CHKERRQ(ierr);
12249e5bc791SBarry Smith   *outits = oits;
12259e5bc791SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
12269e5bc791SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
12279e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
12289e5bc791SBarry Smith   ierr = HYPRE_StructPFMGSetMaxIter(jac->hsolver,jac->its);CHKERRQ(ierr);
12299e5bc791SBarry Smith   PetscFunctionReturn(0);
12309e5bc791SBarry Smith }
12319e5bc791SBarry Smith 
12329e5bc791SBarry Smith 
1233ebc551c0SBarry Smith 
1234ebc551c0SBarry Smith /*MC
1235ebc551c0SBarry Smith      PCPFMG - the hypre PFMG multigrid solver
1236ebc551c0SBarry Smith 
1237ebc551c0SBarry Smith    Level: advanced
1238ebc551c0SBarry Smith 
12399e5bc791SBarry Smith    Options Database:
12409e5bc791SBarry Smith + -pc_pfmg_its <its> number of iterations of PFMG to use as preconditioner
12419e5bc791SBarry Smith . -pc_pfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
12429e5bc791SBarry Smith . -pc_pfmg_num_post_relax <steps> number of smoothing steps after coarse grid
12439e5bc791SBarry Smith . -pc_pfmg_tol <tol> tolerance of PFMG
12449e5bc791SBarry 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
12459e5bc791SBarry Smith - -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin
1246f91d8e95SBarry Smith 
12479e5bc791SBarry Smith    Notes:  This is for CELL-centered descretizations
12489e5bc791SBarry Smith 
12499e5bc791SBarry Smith            This must be used with the MATHYPRESTRUCT matrix type
12509e5bc791SBarry Smith .
12519e5bc791SBarry Smith            This is less general than in hypre, it supports only one block per process defined by a PETSc DA.
12529e5bc791SBarry Smith 
12539e5bc791SBarry Smith .seealso:  PCMG, MATHYPRESTRUCT
1254ebc551c0SBarry Smith M*/
1255ebc551c0SBarry Smith 
1256ebc551c0SBarry Smith EXTERN_C_BEGIN
1257ebc551c0SBarry Smith #undef __FUNCT__
1258ebc551c0SBarry Smith #define __FUNCT__ "PCCreate_PFMG"
1259ebc551c0SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_PFMG(PC pc)
1260ebc551c0SBarry Smith {
1261ebc551c0SBarry Smith   PetscErrorCode ierr;
1262ebc551c0SBarry Smith   PC_PFMG        *ex;
1263ebc551c0SBarry Smith 
1264ebc551c0SBarry Smith   PetscFunctionBegin;
1265ebc551c0SBarry Smith   ierr = PetscNew(PC_PFMG,&ex);CHKERRQ(ierr);\
126668326731SBarry Smith   pc->data = ex;
1267ebc551c0SBarry Smith 
12689e5bc791SBarry Smith   ex->its            = 1;
12699e5bc791SBarry Smith   ex->tol            = 1.e-8;
12709e5bc791SBarry Smith   ex->relax_type     = 1;
12719e5bc791SBarry Smith   ex->rap_type       = 0;
12729e5bc791SBarry Smith   ex->num_pre_relax  = 1;
12739e5bc791SBarry Smith   ex->num_post_relax = 1;
12749e5bc791SBarry Smith 
1275ebc551c0SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_PFMG;
1276ebc551c0SBarry Smith   pc->ops->view            = PCView_PFMG;
1277ebc551c0SBarry Smith   pc->ops->destroy         = PCDestroy_PFMG;
1278f91d8e95SBarry Smith   pc->ops->apply           = PCApply_PFMG;
12799e5bc791SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_PFMG;
128068326731SBarry Smith   pc->ops->setup           = PCSetUp_PFMG;
1281f91d8e95SBarry Smith   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr);
128268326731SBarry Smith   ierr = HYPRE_StructPFMGCreate(ex->hcomm,&ex->hsolver);CHKERRQ(ierr);
1283ebc551c0SBarry Smith   PetscFunctionReturn(0);
1284ebc551c0SBarry Smith }
1285ebc551c0SBarry Smith EXTERN_C_END
1286