xref: /petsc/src/ksp/pc/impls/hypre/hypre.c (revision 390e7148abbd3fdb3c6aa2d3457154909d3a6cea)
116d9e3a6SLisandro Dalcin #define PETSCKSP_DLL
216d9e3a6SLisandro Dalcin 
316d9e3a6SLisandro Dalcin /*
416d9e3a6SLisandro Dalcin    Provides an interface to the LLNL package hypre
516d9e3a6SLisandro Dalcin */
616d9e3a6SLisandro Dalcin #include "private/pcimpl.h"          /*I "petscpc.h" I*/
716d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
816d9e3a6SLisandro Dalcin #include "HYPRE.h"
916d9e3a6SLisandro Dalcin #include "IJ_mv.h"
1016d9e3a6SLisandro Dalcin #include "parcsr_ls.h"
1116d9e3a6SLisandro Dalcin EXTERN_C_END
1216d9e3a6SLisandro Dalcin 
1316d9e3a6SLisandro Dalcin EXTERN PetscErrorCode MatHYPRE_IJMatrixCreate(Mat,HYPRE_IJMatrix*);
1416d9e3a6SLisandro Dalcin EXTERN PetscErrorCode MatHYPRE_IJMatrixCopy(Mat,HYPRE_IJMatrix);
1516d9e3a6SLisandro Dalcin EXTERN PetscErrorCode VecHYPRE_IJVectorCreate(Vec,HYPRE_IJVector*);
1616d9e3a6SLisandro Dalcin 
1716d9e3a6SLisandro Dalcin /*
1816d9e3a6SLisandro Dalcin    Private context (data structure) for the  preconditioner.
1916d9e3a6SLisandro Dalcin */
2016d9e3a6SLisandro Dalcin typedef struct {
2116d9e3a6SLisandro Dalcin   HYPRE_Solver       hsolver;
2216d9e3a6SLisandro Dalcin   HYPRE_IJMatrix     ij;
2316d9e3a6SLisandro Dalcin   HYPRE_IJVector     b,x;
2416d9e3a6SLisandro Dalcin 
2516d9e3a6SLisandro Dalcin   PetscErrorCode     (*destroy)(HYPRE_Solver);
2616d9e3a6SLisandro Dalcin   PetscErrorCode     (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
2716d9e3a6SLisandro Dalcin   PetscErrorCode     (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
2816d9e3a6SLisandro Dalcin 
2916d9e3a6SLisandro Dalcin   MPI_Comm           comm_hypre;
3016d9e3a6SLisandro Dalcin   char              *hypre_type;
3116d9e3a6SLisandro Dalcin 
3216d9e3a6SLisandro Dalcin   /* options for Pilut and BoomerAMG*/
3316d9e3a6SLisandro Dalcin   int                maxiter;
3416d9e3a6SLisandro Dalcin   double             tol;
3516d9e3a6SLisandro Dalcin 
3616d9e3a6SLisandro Dalcin   /* options for Pilut */
3716d9e3a6SLisandro Dalcin   int                factorrowsize;
3816d9e3a6SLisandro Dalcin 
3916d9e3a6SLisandro Dalcin   /* options for ParaSails */
4016d9e3a6SLisandro Dalcin   int                nlevels;
4116d9e3a6SLisandro Dalcin   double             threshhold;
4216d9e3a6SLisandro Dalcin   double             filter;
4316d9e3a6SLisandro Dalcin   int                sym;
4416d9e3a6SLisandro Dalcin   double             loadbal;
4516d9e3a6SLisandro Dalcin   int                logging;
4616d9e3a6SLisandro Dalcin   int                ruse;
4716d9e3a6SLisandro Dalcin   int                symt;
4816d9e3a6SLisandro Dalcin 
4916d9e3a6SLisandro Dalcin   /* options for Euclid */
5016d9e3a6SLisandro Dalcin   PetscTruth         bjilu;
5116d9e3a6SLisandro Dalcin   int                levels;
5216d9e3a6SLisandro Dalcin 
5316d9e3a6SLisandro Dalcin   /* options for Euclid and BoomerAMG */
5416d9e3a6SLisandro Dalcin   PetscTruth         printstatistics;
5516d9e3a6SLisandro Dalcin 
5616d9e3a6SLisandro Dalcin   /* options for BoomerAMG */
5716d9e3a6SLisandro Dalcin   int                cycletype;
5816d9e3a6SLisandro Dalcin   int                maxlevels;
5916d9e3a6SLisandro Dalcin   double             strongthreshold;
6016d9e3a6SLisandro Dalcin   double             maxrowsum;
6116d9e3a6SLisandro Dalcin   int                gridsweeps[4];
6216d9e3a6SLisandro Dalcin   int                coarsentype;
6316d9e3a6SLisandro Dalcin   int                measuretype;
6416d9e3a6SLisandro Dalcin   int                relaxtype[4];
6516d9e3a6SLisandro Dalcin   double             relaxweight;
6616d9e3a6SLisandro Dalcin   double             outerrelaxweight;
6716d9e3a6SLisandro Dalcin   int                relaxorder;
6816d9e3a6SLisandro Dalcin   double             truncfactor;
6916d9e3a6SLisandro Dalcin   PetscTruth         applyrichardson;
7016d9e3a6SLisandro Dalcin 
7116d9e3a6SLisandro Dalcin } PC_HYPRE;
7216d9e3a6SLisandro Dalcin 
7316d9e3a6SLisandro Dalcin 
7416d9e3a6SLisandro Dalcin #undef __FUNCT__
7516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetUp_HYPRE"
7616d9e3a6SLisandro Dalcin static PetscErrorCode PCSetUp_HYPRE(PC pc)
7716d9e3a6SLisandro Dalcin {
7816d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
7916d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
8016d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
8116d9e3a6SLisandro Dalcin   HYPRE_ParVector    bv,xv;
8216d9e3a6SLisandro Dalcin   PetscInt           bs;
8316d9e3a6SLisandro Dalcin   int                hierr;
8416d9e3a6SLisandro Dalcin 
8516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
8616d9e3a6SLisandro Dalcin   if (!jac->hypre_type) {
8716d9e3a6SLisandro Dalcin     ierr = PCHYPRESetType(pc,"pilut");CHKERRQ(ierr);
8816d9e3a6SLisandro Dalcin   }
8916d9e3a6SLisandro Dalcin #if 1
9016d9e3a6SLisandro Dalcin   if (!pc->setupcalled) {
9116d9e3a6SLisandro Dalcin     /* create the matrix and vectors the first time through */
9216d9e3a6SLisandro Dalcin     Vec x,b;
9316d9e3a6SLisandro Dalcin     ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr);
9416d9e3a6SLisandro Dalcin     ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr);
9516d9e3a6SLisandro Dalcin     ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr);
9616d9e3a6SLisandro Dalcin     ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr);
9716d9e3a6SLisandro Dalcin     ierr = VecDestroy(x);CHKERRQ(ierr);
9816d9e3a6SLisandro Dalcin     ierr = VecDestroy(b);CHKERRQ(ierr);
9916d9e3a6SLisandro Dalcin   } else if (pc->flag != SAME_NONZERO_PATTERN) {
10016d9e3a6SLisandro Dalcin     /* rebuild the matrix from scratch */
10116d9e3a6SLisandro Dalcin     ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr);
10216d9e3a6SLisandro Dalcin     ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr);
10316d9e3a6SLisandro Dalcin   }
10416d9e3a6SLisandro Dalcin #else
10516d9e3a6SLisandro Dalcin   if (!jac->ij) { /* create the matrix the first time through */
10616d9e3a6SLisandro Dalcin     ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr);
10716d9e3a6SLisandro Dalcin   }
10816d9e3a6SLisandro Dalcin   if (!jac->b) { /* create the vectors the first time through */
10916d9e3a6SLisandro Dalcin     Vec x,b;
11016d9e3a6SLisandro Dalcin     ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr);
11116d9e3a6SLisandro Dalcin     ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr);
11216d9e3a6SLisandro Dalcin     ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr);
11316d9e3a6SLisandro Dalcin     ierr = VecDestroy(x);CHKERRQ(ierr);
11416d9e3a6SLisandro Dalcin     ierr = VecDestroy(b);CHKERRQ(ierr);
11516d9e3a6SLisandro Dalcin   }
11616d9e3a6SLisandro Dalcin #endif
11716d9e3a6SLisandro Dalcin   /* special case for BoomerAMG */
11816d9e3a6SLisandro Dalcin   if (jac->setup == HYPRE_BoomerAMGSetup) {
11916d9e3a6SLisandro Dalcin     ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr);
12016d9e3a6SLisandro Dalcin     if (bs > 1) {
12116d9e3a6SLisandro Dalcin       ierr = HYPRE_BoomerAMGSetNumFunctions(jac->hsolver,bs);CHKERRQ(ierr);
12216d9e3a6SLisandro Dalcin     }
12316d9e3a6SLisandro Dalcin   };
12416d9e3a6SLisandro Dalcin   ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr);
12516d9e3a6SLisandro Dalcin   ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr);
12616d9e3a6SLisandro Dalcin   ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&bv);CHKERRQ(ierr);
12716d9e3a6SLisandro Dalcin   ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&xv);CHKERRQ(ierr);
12816d9e3a6SLisandro Dalcin   hierr = (*jac->setup)(jac->hsolver,hmat,bv,xv);
12916d9e3a6SLisandro Dalcin   if (hierr) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE setup, error code %d",hierr);
13016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
13116d9e3a6SLisandro Dalcin }
13216d9e3a6SLisandro Dalcin 
13316d9e3a6SLisandro Dalcin /*
13416d9e3a6SLisandro Dalcin     Replaces the address where the HYPRE vector points to its data with the address of
13516d9e3a6SLisandro Dalcin   PETSc's data. Saves the old address so it can be reset when we are finished with it.
13616d9e3a6SLisandro Dalcin   Allows use to get the data into a HYPRE vector without the cost of memcopies
13716d9e3a6SLisandro Dalcin */
13816d9e3a6SLisandro Dalcin #define HYPREReplacePointer(b,newvalue,savedvalue) {\
13916d9e3a6SLisandro Dalcin    hypre_ParVector *par_vector   = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));\
14016d9e3a6SLisandro Dalcin    hypre_Vector    *local_vector = hypre_ParVectorLocalVector(par_vector);\
14116d9e3a6SLisandro Dalcin    savedvalue         = local_vector->data;\
14216d9e3a6SLisandro Dalcin    local_vector->data = newvalue;}
14316d9e3a6SLisandro Dalcin 
14416d9e3a6SLisandro Dalcin #undef __FUNCT__
14516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApply_HYPRE"
14616d9e3a6SLisandro Dalcin static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
14716d9e3a6SLisandro Dalcin {
14816d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
14916d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
15016d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
15116d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
15216d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
15316d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
15416d9e3a6SLisandro Dalcin   int                hierr;
15516d9e3a6SLisandro Dalcin 
15616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
15716d9e3a6SLisandro Dalcin   if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);}
15816d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
15916d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
16016d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
16116d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
16216d9e3a6SLisandro Dalcin 
16316d9e3a6SLisandro Dalcin   ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr);
16416d9e3a6SLisandro Dalcin   ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr);
16516d9e3a6SLisandro Dalcin   ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr);
16616d9e3a6SLisandro Dalcin   hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv);
16716d9e3a6SLisandro Dalcin   /* error code of 1 in BoomerAMG merely means convergence not achieved */
16816d9e3a6SLisandro Dalcin   if (hierr && (hierr != 1 || jac->solve != HYPRE_BoomerAMGSolve)) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
16916d9e3a6SLisandro Dalcin   if (hierr) hypre__global_error = 0;
17016d9e3a6SLisandro Dalcin 
17116d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
17216d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
17316d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
17416d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
17516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
17616d9e3a6SLisandro Dalcin }
17716d9e3a6SLisandro Dalcin 
17816d9e3a6SLisandro Dalcin #undef __FUNCT__
17916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCDestroy_HYPRE"
18016d9e3a6SLisandro Dalcin static PetscErrorCode PCDestroy_HYPRE(PC pc)
18116d9e3a6SLisandro Dalcin {
18216d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
18316d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
18416d9e3a6SLisandro Dalcin 
18516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
18616d9e3a6SLisandro Dalcin   if (jac->ij) { ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); }
18716d9e3a6SLisandro Dalcin   if (jac->b)  { ierr = HYPRE_IJVectorDestroy(jac->b);CHKERRQ(ierr);  }
18816d9e3a6SLisandro Dalcin   if (jac->x)  { ierr = HYPRE_IJVectorDestroy(jac->x);CHKERRQ(ierr);  }
18916d9e3a6SLisandro Dalcin   if (jac->destroy) { ierr = (*jac->destroy)(jac->hsolver);CHKERRQ(ierr); }
19016d9e3a6SLisandro Dalcin   ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr);
19116d9e3a6SLisandro Dalcin   if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);}
19216d9e3a6SLisandro Dalcin   ierr = PetscFree(jac);CHKERRQ(ierr);
19316d9e3a6SLisandro Dalcin 
19416d9e3a6SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr);
19516d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);CHKERRQ(ierr);
19616d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);CHKERRQ(ierr);
19716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
19816d9e3a6SLisandro Dalcin }
19916d9e3a6SLisandro Dalcin 
20016d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
20116d9e3a6SLisandro Dalcin #undef __FUNCT__
20216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut"
20316d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc)
20416d9e3a6SLisandro Dalcin {
20516d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
20616d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
20716d9e3a6SLisandro Dalcin   PetscTruth     flag;
20816d9e3a6SLisandro Dalcin 
20916d9e3a6SLisandro Dalcin   PetscFunctionBegin;
21016d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr);
21116d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr);
21216d9e3a6SLisandro Dalcin   if (flag) {
21316d9e3a6SLisandro Dalcin     ierr = HYPRE_ParCSRPilutSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
21416d9e3a6SLisandro Dalcin   }
21516d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr);
21616d9e3a6SLisandro Dalcin   if (flag) {
21716d9e3a6SLisandro Dalcin     ierr = HYPRE_ParCSRPilutSetDropTolerance(jac->hsolver,jac->tol);CHKERRQ(ierr);
21816d9e3a6SLisandro Dalcin   }
21916d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr);
22016d9e3a6SLisandro Dalcin   if (flag) {
22116d9e3a6SLisandro Dalcin     ierr = HYPRE_ParCSRPilutSetFactorRowSize(jac->hsolver,jac->factorrowsize);CHKERRQ(ierr);
22216d9e3a6SLisandro Dalcin   }
22316d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
22416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
22516d9e3a6SLisandro Dalcin }
22616d9e3a6SLisandro Dalcin 
22716d9e3a6SLisandro Dalcin #undef __FUNCT__
22816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Pilut"
22916d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
23016d9e3a6SLisandro Dalcin {
23116d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
23216d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
23316d9e3a6SLisandro Dalcin   PetscTruth     iascii;
23416d9e3a6SLisandro Dalcin 
23516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
23616d9e3a6SLisandro Dalcin   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
23716d9e3a6SLisandro Dalcin   if (iascii) {
23816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");CHKERRQ(ierr);
23916d9e3a6SLisandro Dalcin     if (jac->maxiter != PETSC_DEFAULT) {
24016d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr);
24116d9e3a6SLisandro Dalcin     } else {
24216d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr);
24316d9e3a6SLisandro Dalcin     }
24416d9e3a6SLisandro Dalcin     if (jac->tol != PETSC_DEFAULT) {
24516d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr);
24616d9e3a6SLisandro Dalcin     } else {
24716d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr);
24816d9e3a6SLisandro Dalcin     }
24916d9e3a6SLisandro Dalcin     if (jac->factorrowsize != PETSC_DEFAULT) {
25016d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr);
25116d9e3a6SLisandro Dalcin     } else {
25216d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");CHKERRQ(ierr);
25316d9e3a6SLisandro Dalcin     }
25416d9e3a6SLisandro Dalcin   }
25516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
25616d9e3a6SLisandro Dalcin }
25716d9e3a6SLisandro Dalcin 
25816d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
25916d9e3a6SLisandro Dalcin #undef __FUNCT__
26016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid"
26116d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc)
26216d9e3a6SLisandro Dalcin {
26316d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
26416d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
26516d9e3a6SLisandro Dalcin   PetscTruth     flag;
266*390e7148SBarry Smith   char           *args[8],levels[16];
267*390e7148SBarry Smith   PetscInt       cnt = 0;
26816d9e3a6SLisandro Dalcin 
26916d9e3a6SLisandro Dalcin   PetscFunctionBegin;
27016d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr);
27116d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr);
27216d9e3a6SLisandro Dalcin   if (flag) {
27316d9e3a6SLisandro Dalcin     if (jac->levels < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels);
27416d9e3a6SLisandro Dalcin     sprintf(levels,"%d",jac->levels);
275*390e7148SBarry Smith     args[cnt++] = (char*)"-level"; args[cnt++] = levels;
27616d9e3a6SLisandro Dalcin   }
27716d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);CHKERRQ(ierr);
27816d9e3a6SLisandro Dalcin   if (jac->bjilu) {
279*390e7148SBarry Smith     args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1";
28016d9e3a6SLisandro Dalcin   }
28116d9e3a6SLisandro Dalcin 
28216d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);CHKERRQ(ierr);
28316d9e3a6SLisandro Dalcin   if (jac->printstatistics) {
284*390e7148SBarry Smith     args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1";
285*390e7148SBarry Smith     args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1";
28616d9e3a6SLisandro Dalcin   }
28716d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
288*390e7148SBarry Smith   if (cnt) {
289*390e7148SBarry Smith     ierr = HYPRE_EuclidSetParams(jac->hsolver,cnt,args);CHKERRQ(ierr);
290*390e7148SBarry Smith   }
29116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
29216d9e3a6SLisandro Dalcin }
29316d9e3a6SLisandro Dalcin 
29416d9e3a6SLisandro Dalcin #undef __FUNCT__
29516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Euclid"
29616d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
29716d9e3a6SLisandro Dalcin {
29816d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
29916d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
30016d9e3a6SLisandro Dalcin   PetscTruth     iascii;
30116d9e3a6SLisandro Dalcin 
30216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
30316d9e3a6SLisandro Dalcin   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
30416d9e3a6SLisandro Dalcin   if (iascii) {
30516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n");CHKERRQ(ierr);
30616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr);
30716d9e3a6SLisandro Dalcin     if (jac->bjilu) {
30816d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr);
30916d9e3a6SLisandro Dalcin     }
31016d9e3a6SLisandro Dalcin   }
31116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
31216d9e3a6SLisandro Dalcin }
31316d9e3a6SLisandro Dalcin 
31416d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
31516d9e3a6SLisandro Dalcin 
31616d9e3a6SLisandro Dalcin #undef __FUNCT__
31716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG"
31816d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
31916d9e3a6SLisandro Dalcin {
32016d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
32116d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
32216d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
32316d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
32416d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
32516d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
32616d9e3a6SLisandro Dalcin   int                hierr;
32716d9e3a6SLisandro Dalcin 
32816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
32916d9e3a6SLisandro Dalcin   ierr = VecSet(x,0.0);CHKERRQ(ierr);
33016d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
33116d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
33216d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
33316d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
33416d9e3a6SLisandro Dalcin 
33516d9e3a6SLisandro Dalcin   ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr);
33616d9e3a6SLisandro Dalcin   ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr);
33716d9e3a6SLisandro Dalcin   ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr);
33816d9e3a6SLisandro Dalcin 
33916d9e3a6SLisandro Dalcin   hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv);
34016d9e3a6SLisandro Dalcin   /* error code of 1 in BoomerAMG merely means convergence not achieved */
34116d9e3a6SLisandro Dalcin   if (hierr && (hierr != 1)) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
34216d9e3a6SLisandro Dalcin   if (hierr) hypre__global_error = 0;
34316d9e3a6SLisandro Dalcin 
34416d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
34516d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
34616d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
34716d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
34816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
34916d9e3a6SLisandro Dalcin }
35016d9e3a6SLisandro Dalcin 
35116d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
35216d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout"};
35316d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
35416d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","","SOR/Jacobi","backward-SOR/Jacobi","","symmetric-SOR/Jacobi","","","Gaussian-elimination"};
35516d9e3a6SLisandro Dalcin 
35616d9e3a6SLisandro Dalcin #undef __FUNCT__
35716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG"
35816d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
35916d9e3a6SLisandro Dalcin {
36016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
36116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
36216d9e3a6SLisandro Dalcin   int            n,indx;
36316d9e3a6SLisandro Dalcin   PetscTruth     flg, tmp_truth;
36416d9e3a6SLisandro Dalcin   double         tmpdbl, twodbl[2];
36516d9e3a6SLisandro Dalcin 
36616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
36716d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr);
36816d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr);
36916d9e3a6SLisandro Dalcin   if (flg) {
37016d9e3a6SLisandro Dalcin     jac->cycletype = indx;
37116d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr);
37216d9e3a6SLisandro Dalcin   }
37316d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr);
37416d9e3a6SLisandro Dalcin   if (flg) {
37516d9e3a6SLisandro Dalcin     if (jac->maxlevels < 2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
37616d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr);
37716d9e3a6SLisandro Dalcin   }
37816d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr);
37916d9e3a6SLisandro Dalcin   if (flg) {
38016d9e3a6SLisandro Dalcin     if (jac->maxiter < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
38116d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
38216d9e3a6SLisandro Dalcin   }
38316d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_tol","Convergence tolerance PER hypre call","None",jac->tol,&jac->tol,&flg);CHKERRQ(ierr);
38416d9e3a6SLisandro Dalcin   if (flg) {
38516d9e3a6SLisandro Dalcin     if (jac->tol < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be great than or equal zero",jac->tol);
38616d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
38716d9e3a6SLisandro Dalcin   }
38816d9e3a6SLisandro Dalcin 
38916d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr);
39016d9e3a6SLisandro Dalcin   if (flg) {
39116d9e3a6SLisandro Dalcin     if (jac->truncfactor < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor);
39216d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr);
39316d9e3a6SLisandro Dalcin   }
39416d9e3a6SLisandro Dalcin 
39516d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr);
39616d9e3a6SLisandro Dalcin   if (flg) {
39716d9e3a6SLisandro Dalcin     if (jac->strongthreshold < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold);
39816d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr);
39916d9e3a6SLisandro Dalcin   }
40016d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr);
40116d9e3a6SLisandro Dalcin   if (flg) {
40216d9e3a6SLisandro Dalcin     if (jac->maxrowsum < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum);
40316d9e3a6SLisandro Dalcin     if (jac->maxrowsum > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum);
40416d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr);
40516d9e3a6SLisandro Dalcin   }
40616d9e3a6SLisandro Dalcin 
40716d9e3a6SLisandro Dalcin   /* Grid sweeps */
40816d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_all","Number of sweeps for all grid levels (fine, up, and down)","None",jac->gridsweeps[0],&indx,&flg);CHKERRQ(ierr);
40916d9e3a6SLisandro Dalcin   if (flg) {
41016d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver,indx);CHKERRQ(ierr);
41116d9e3a6SLisandro Dalcin     /* modify the jac structure so we can view the updated options with PC_View */
41216d9e3a6SLisandro Dalcin     jac->gridsweeps[0] = indx;
41316d9e3a6SLisandro Dalcin     jac->gridsweeps[1] = jac->gridsweeps[2] = jac->gridsweeps[0];
41416d9e3a6SLisandro Dalcin     jac->gridsweeps[3] = 1;  /*The coarse level is not affected by this function - hypre code sets to 1*/
41516d9e3a6SLisandro Dalcin   }
41616d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_fine","Number of sweeps for the fine level","None",jac->gridsweeps[0],&indx,&flg);CHKERRQ(ierr);
41716d9e3a6SLisandro Dalcin   if (flg) {
41816d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 0);CHKERRQ(ierr);
41916d9e3a6SLisandro Dalcin     jac->gridsweeps[0] = indx;
42016d9e3a6SLisandro Dalcin   }
42116d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr);
42216d9e3a6SLisandro Dalcin   if (flg) {
42316d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 1);CHKERRQ(ierr);
42416d9e3a6SLisandro Dalcin     jac->gridsweeps[1] = indx;
42516d9e3a6SLisandro Dalcin   }
42616d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr);
42716d9e3a6SLisandro Dalcin   if (flg) {
42816d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 2);CHKERRQ(ierr);
42916d9e3a6SLisandro Dalcin     jac->gridsweeps[2] = indx;
43016d9e3a6SLisandro Dalcin   }
43116d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[3],&indx,&flg);CHKERRQ(ierr);
43216d9e3a6SLisandro Dalcin   if (flg) {
43316d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 3);CHKERRQ(ierr);
43416d9e3a6SLisandro Dalcin     jac->gridsweeps[3] = indx;
43516d9e3a6SLisandro Dalcin   }
43616d9e3a6SLisandro Dalcin 
43716d9e3a6SLisandro Dalcin   /* Relax type */
43816d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_all","Relax type for fine, up, and down cycles (coarse level set to gaussian elimination)","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);CHKERRQ(ierr);
43916d9e3a6SLisandro Dalcin   if (flg) {
44016d9e3a6SLisandro Dalcin     jac->relaxtype[0] = jac->relaxtype[1] = jac->relaxtype[2] = indx;
44116d9e3a6SLisandro Dalcin     jac->relaxtype[3] = 9; /* hypre code sets coarse grid to 9 (G.E.)*/
44216d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, indx);CHKERRQ(ierr);
44316d9e3a6SLisandro Dalcin   }
44416d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_fine","Relax type on fine grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);CHKERRQ(ierr);
44516d9e3a6SLisandro Dalcin   if (flg) {
44616d9e3a6SLisandro Dalcin     jac->relaxtype[0] = indx;
44716d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 0);CHKERRQ(ierr);
44816d9e3a6SLisandro Dalcin   }
44916d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);CHKERRQ(ierr);
45016d9e3a6SLisandro Dalcin   if (flg) {
45116d9e3a6SLisandro Dalcin     jac->relaxtype[1] = indx;
45216d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 1);CHKERRQ(ierr);
45316d9e3a6SLisandro Dalcin   }
45416d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);CHKERRQ(ierr);
45516d9e3a6SLisandro Dalcin   if (flg) {
45616d9e3a6SLisandro Dalcin     jac->relaxtype[2] = indx;
45716d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 2);CHKERRQ(ierr);
45816d9e3a6SLisandro Dalcin   }
45916d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr);
46016d9e3a6SLisandro Dalcin   if (flg) {
46116d9e3a6SLisandro Dalcin     jac->relaxtype[3] = indx;
46216d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 3);CHKERRQ(ierr);
46316d9e3a6SLisandro Dalcin   }
46416d9e3a6SLisandro Dalcin 
46516d9e3a6SLisandro Dalcin   /* Relaxation Weight */
46616d9e3a6SLisandro 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);
46716d9e3a6SLisandro Dalcin   if (flg) {
46816d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetRelaxWt(jac->hsolver,tmpdbl);CHKERRQ(ierr);
46916d9e3a6SLisandro Dalcin     jac->relaxweight = tmpdbl;
47016d9e3a6SLisandro Dalcin   }
47116d9e3a6SLisandro Dalcin 
47216d9e3a6SLisandro Dalcin   n=2;
47316d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
47416d9e3a6SLisandro 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);
47516d9e3a6SLisandro Dalcin   if (flg) {
47616d9e3a6SLisandro Dalcin     if (n == 2) {
47716d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
47816d9e3a6SLisandro Dalcin       ierr = HYPRE_BoomerAMGSetLevelRelaxWt(jac->hsolver,twodbl[0],indx);CHKERRQ(ierr);
47916d9e3a6SLisandro Dalcin     } else {
48016d9e3a6SLisandro 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);
48116d9e3a6SLisandro Dalcin     }
48216d9e3a6SLisandro Dalcin   }
48316d9e3a6SLisandro Dalcin 
48416d9e3a6SLisandro Dalcin   /* Outer relaxation Weight */
48516d9e3a6SLisandro 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);
48616d9e3a6SLisandro Dalcin   if (flg) {
48716d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetOuterWt( jac->hsolver, tmpdbl);CHKERRQ(ierr);
48816d9e3a6SLisandro Dalcin     jac->outerrelaxweight = tmpdbl;
48916d9e3a6SLisandro Dalcin   }
49016d9e3a6SLisandro Dalcin 
49116d9e3a6SLisandro Dalcin   n=2;
49216d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
49316d9e3a6SLisandro 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);
49416d9e3a6SLisandro Dalcin   if (flg) {
49516d9e3a6SLisandro Dalcin     if (n == 2) {
49616d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
49716d9e3a6SLisandro Dalcin       ierr = HYPRE_BoomerAMGSetLevelOuterWt( jac->hsolver, twodbl[0], indx);CHKERRQ(ierr);
49816d9e3a6SLisandro Dalcin     } else {
49916d9e3a6SLisandro 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);
50016d9e3a6SLisandro Dalcin     }
50116d9e3a6SLisandro Dalcin   }
50216d9e3a6SLisandro Dalcin 
50316d9e3a6SLisandro Dalcin   /* the Relax Order */
50416d9e3a6SLisandro Dalcin   /* ierr = PetscOptionsName("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", &flg);CHKERRQ(ierr); */
50516d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
50616d9e3a6SLisandro Dalcin 
50716d9e3a6SLisandro Dalcin   if (flg) {
50816d9e3a6SLisandro Dalcin     jac->relaxorder = 0;
50916d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr);
51016d9e3a6SLisandro Dalcin   }
51116d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr);
51216d9e3a6SLisandro Dalcin   if (flg) {
51316d9e3a6SLisandro Dalcin     jac->measuretype = indx;
51416d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr);
51516d9e3a6SLisandro Dalcin   }
51616d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,7,HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr);
51716d9e3a6SLisandro Dalcin   if (flg) {
51816d9e3a6SLisandro Dalcin     jac->coarsentype = indx;
51916d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr);
52016d9e3a6SLisandro Dalcin   }
52116d9e3a6SLisandro Dalcin   ierr = PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);CHKERRQ(ierr);
52216d9e3a6SLisandro Dalcin   if (flg) {
52316d9e3a6SLisandro Dalcin     int level=3;
52416d9e3a6SLisandro Dalcin     jac->printstatistics = PETSC_TRUE;
52516d9e3a6SLisandro Dalcin     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
52616d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetPrintLevel(jac->hsolver,level);CHKERRQ(ierr);
52716d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetDebugFlag(jac->hsolver,level);CHKERRQ(ierr);
52816d9e3a6SLisandro Dalcin   }
52916d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
53016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
53116d9e3a6SLisandro Dalcin }
53216d9e3a6SLisandro Dalcin 
53316d9e3a6SLisandro Dalcin #undef __FUNCT__
53416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG"
53516d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its)
53616d9e3a6SLisandro Dalcin {
53716d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
53816d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
53916d9e3a6SLisandro Dalcin 
54016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
54116d9e3a6SLisandro Dalcin   ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,its*jac->maxiter);CHKERRQ(ierr);
54216d9e3a6SLisandro Dalcin   ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,PetscMin(rtol,jac->tol));CHKERRQ(ierr);
54316d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_TRUE;
54416d9e3a6SLisandro Dalcin   ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr);
54516d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_FALSE;
54616d9e3a6SLisandro Dalcin   ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
54716d9e3a6SLisandro Dalcin   ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
54816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
54916d9e3a6SLisandro Dalcin }
55016d9e3a6SLisandro Dalcin 
55116d9e3a6SLisandro Dalcin 
55216d9e3a6SLisandro Dalcin #undef __FUNCT__
55316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_BoomerAMG"
55416d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
55516d9e3a6SLisandro Dalcin {
55616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
55716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
55816d9e3a6SLisandro Dalcin   PetscTruth     iascii;
55916d9e3a6SLisandro Dalcin 
56016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
56116d9e3a6SLisandro Dalcin   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
56216d9e3a6SLisandro Dalcin   if (iascii) {
56316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr);
56416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr);
56516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr);
56616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr);
56716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr);
56816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr);
56916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr);
57016d9e3a6SLisandro Dalcin 
57116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on fine grid %d\n",jac->gridsweeps[0]);CHKERRQ(ierr);
57216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps down         %d\n",jac->gridsweeps[1]);CHKERRQ(ierr);
57316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps up           %d\n",jac->gridsweeps[2]);CHKERRQ(ierr);
57416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on coarse    %d\n",jac->gridsweeps[3]);CHKERRQ(ierr);
57516d9e3a6SLisandro Dalcin 
57616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on fine grid  %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr);
57716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr);
57816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr);
57916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[3]]);CHKERRQ(ierr);
58016d9e3a6SLisandro Dalcin 
58116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax weight  (all)      %G\n",jac->relaxweight);CHKERRQ(ierr);
58216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr);
58316d9e3a6SLisandro Dalcin 
58416d9e3a6SLisandro Dalcin     if (jac->relaxorder) {
58516d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr);
58616d9e3a6SLisandro Dalcin     } else {
58716d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr);
58816d9e3a6SLisandro Dalcin     }
58916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr);
59016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr);
59116d9e3a6SLisandro Dalcin   }
59216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
59316d9e3a6SLisandro Dalcin }
59416d9e3a6SLisandro Dalcin 
59516d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
59616d9e3a6SLisandro Dalcin #undef __FUNCT__
59716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails"
59816d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc)
59916d9e3a6SLisandro Dalcin {
60016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
60116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
60216d9e3a6SLisandro Dalcin   int            indx;
60316d9e3a6SLisandro Dalcin   PetscTruth     flag;
60416d9e3a6SLisandro Dalcin   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};
60516d9e3a6SLisandro Dalcin 
60616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
60716d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr);
60816d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr);
60916d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr);
61016d9e3a6SLisandro Dalcin   if (flag) {
61116d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr);
61216d9e3a6SLisandro Dalcin   }
61316d9e3a6SLisandro Dalcin 
61416d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr);
61516d9e3a6SLisandro Dalcin   if (flag) {
61616d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr);
61716d9e3a6SLisandro Dalcin   }
61816d9e3a6SLisandro Dalcin 
61916d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr);
62016d9e3a6SLisandro Dalcin   if (flag) {
62116d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr);
62216d9e3a6SLisandro Dalcin   }
62316d9e3a6SLisandro Dalcin 
62416d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,&flag);CHKERRQ(ierr);
62516d9e3a6SLisandro Dalcin   if (flag) {
62616d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr);
62716d9e3a6SLisandro Dalcin   }
62816d9e3a6SLisandro Dalcin 
62916d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,&flag);CHKERRQ(ierr);
63016d9e3a6SLisandro Dalcin   if (flag) {
63116d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr);
63216d9e3a6SLisandro Dalcin   }
63316d9e3a6SLisandro Dalcin 
63416d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);CHKERRQ(ierr);
63516d9e3a6SLisandro Dalcin   if (flag) {
63616d9e3a6SLisandro Dalcin     jac->symt = indx;
63716d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr);
63816d9e3a6SLisandro Dalcin   }
63916d9e3a6SLisandro Dalcin 
64016d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
64116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
64216d9e3a6SLisandro Dalcin }
64316d9e3a6SLisandro Dalcin 
64416d9e3a6SLisandro Dalcin #undef __FUNCT__
64516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_ParaSails"
64616d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
64716d9e3a6SLisandro Dalcin {
64816d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
64916d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
65016d9e3a6SLisandro Dalcin   PetscTruth     iascii;
65116d9e3a6SLisandro Dalcin   const char     *symt = 0;;
65216d9e3a6SLisandro Dalcin 
65316d9e3a6SLisandro Dalcin   PetscFunctionBegin;
65416d9e3a6SLisandro Dalcin   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
65516d9e3a6SLisandro Dalcin   if (iascii) {
65616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");CHKERRQ(ierr);
65716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr);
65816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr);
65916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr);
66016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr);
66116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscTruths[jac->ruse]);CHKERRQ(ierr);
66216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscTruths[jac->logging]);CHKERRQ(ierr);
66316d9e3a6SLisandro Dalcin     if (!jac->symt) {
66416d9e3a6SLisandro Dalcin       symt = "nonsymmetric matrix and preconditioner";
66516d9e3a6SLisandro Dalcin     } else if (jac->symt == 1) {
66616d9e3a6SLisandro Dalcin       symt = "SPD matrix and preconditioner";
66716d9e3a6SLisandro Dalcin     } else if (jac->symt == 2) {
66816d9e3a6SLisandro Dalcin       symt = "nonsymmetric matrix but SPD preconditioner";
66916d9e3a6SLisandro Dalcin     } else {
67016d9e3a6SLisandro Dalcin       SETERRQ1(PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
67116d9e3a6SLisandro Dalcin     }
67216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr);
67316d9e3a6SLisandro Dalcin   }
67416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
67516d9e3a6SLisandro Dalcin }
67616d9e3a6SLisandro Dalcin /* ---------------------------------------------------------------------------------*/
67716d9e3a6SLisandro Dalcin 
67816d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
67916d9e3a6SLisandro Dalcin #undef __FUNCT__
68016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType_HYPRE"
68116d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType_HYPRE(PC pc,const char *name[])
68216d9e3a6SLisandro Dalcin {
68316d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
68416d9e3a6SLisandro Dalcin 
68516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
68616d9e3a6SLisandro Dalcin   *name = jac->hypre_type;
68716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
68816d9e3a6SLisandro Dalcin }
68916d9e3a6SLisandro Dalcin EXTERN_C_END
69016d9e3a6SLisandro Dalcin 
69116d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
69216d9e3a6SLisandro Dalcin #undef __FUNCT__
69316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType_HYPRE"
69416d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType_HYPRE(PC pc,const char name[])
69516d9e3a6SLisandro Dalcin {
69616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
69716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
69816d9e3a6SLisandro Dalcin   PetscTruth     flag;
69916d9e3a6SLisandro Dalcin 
70016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
70116d9e3a6SLisandro Dalcin   if (jac->hypre_type) {
70216d9e3a6SLisandro Dalcin     ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr);
70316d9e3a6SLisandro Dalcin     if (!flag) SETERRQ(PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
70416d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
70516d9e3a6SLisandro Dalcin   } else {
70616d9e3a6SLisandro Dalcin     ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr);
70716d9e3a6SLisandro Dalcin   }
70816d9e3a6SLisandro Dalcin 
70916d9e3a6SLisandro Dalcin   jac->maxiter            = PETSC_DEFAULT;
71016d9e3a6SLisandro Dalcin   jac->tol                = PETSC_DEFAULT;
71116d9e3a6SLisandro Dalcin   jac->printstatistics    = PetscLogPrintInfo;
71216d9e3a6SLisandro Dalcin 
71316d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr);
71416d9e3a6SLisandro Dalcin   if (flag) {
71516d9e3a6SLisandro Dalcin     ierr                    = HYPRE_ParCSRPilutCreate(jac->comm_hypre,&jac->hsolver);
71616d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
71716d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Pilut;
71816d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParCSRPilutDestroy;
71916d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParCSRPilutSetup;
72016d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParCSRPilutSolve;
72116d9e3a6SLisandro Dalcin     jac->factorrowsize      = PETSC_DEFAULT;
72216d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
72316d9e3a6SLisandro Dalcin   }
72416d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr);
72516d9e3a6SLisandro Dalcin   if (flag) {
72616d9e3a6SLisandro Dalcin     ierr                    = HYPRE_ParaSailsCreate(jac->comm_hypre,&jac->hsolver);
72716d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
72816d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_ParaSails;
72916d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParaSailsDestroy;
73016d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParaSailsSetup;
73116d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParaSailsSolve;
73216d9e3a6SLisandro Dalcin     /* initialize */
73316d9e3a6SLisandro Dalcin     jac->nlevels     = 1;
73416d9e3a6SLisandro Dalcin     jac->threshhold  = .1;
73516d9e3a6SLisandro Dalcin     jac->filter      = .1;
73616d9e3a6SLisandro Dalcin     jac->loadbal     = 0;
73716d9e3a6SLisandro Dalcin     if (PetscLogPrintInfo) {
73816d9e3a6SLisandro Dalcin       jac->logging     = (int) PETSC_TRUE;
73916d9e3a6SLisandro Dalcin     } else {
74016d9e3a6SLisandro Dalcin       jac->logging     = (int) PETSC_FALSE;
74116d9e3a6SLisandro Dalcin     }
74216d9e3a6SLisandro Dalcin     jac->ruse = (int) PETSC_FALSE;
74316d9e3a6SLisandro Dalcin     jac->symt   = 0;
74416d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr);
74516d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr);
74616d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr);
74716d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr);
74816d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr);
74916d9e3a6SLisandro Dalcin     ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr);
75016d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
75116d9e3a6SLisandro Dalcin   }
75216d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr);
75316d9e3a6SLisandro Dalcin   if (flag) {
75416d9e3a6SLisandro Dalcin     ierr                    = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver);
75516d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
75616d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Euclid;
75716d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_EuclidDestroy;
75816d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_EuclidSetup;
75916d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_EuclidSolve;
76016d9e3a6SLisandro Dalcin     /* initialization */
76116d9e3a6SLisandro Dalcin     jac->bjilu              = PETSC_FALSE;
76216d9e3a6SLisandro Dalcin     jac->levels             = 1;
76316d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
76416d9e3a6SLisandro Dalcin   }
76516d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr);
76616d9e3a6SLisandro Dalcin   if (flag) {
76716d9e3a6SLisandro Dalcin     ierr                     = HYPRE_BoomerAMGCreate(&jac->hsolver);
76816d9e3a6SLisandro Dalcin     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
76916d9e3a6SLisandro Dalcin     pc->ops->view            = PCView_HYPRE_BoomerAMG;
77016d9e3a6SLisandro Dalcin     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
77116d9e3a6SLisandro Dalcin     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
77216d9e3a6SLisandro Dalcin     jac->destroy             = HYPRE_BoomerAMGDestroy;
77316d9e3a6SLisandro Dalcin     jac->setup               = HYPRE_BoomerAMGSetup;
77416d9e3a6SLisandro Dalcin     jac->solve               = HYPRE_BoomerAMGSolve;
77516d9e3a6SLisandro Dalcin     jac->applyrichardson     = PETSC_FALSE;
77616d9e3a6SLisandro Dalcin     /* these defaults match the hypre defaults */
77716d9e3a6SLisandro Dalcin     jac->cycletype        = 1;
77816d9e3a6SLisandro Dalcin     jac->maxlevels        = 25;
77916d9e3a6SLisandro Dalcin     jac->maxiter          = 1;
78016d9e3a6SLisandro Dalcin     jac->tol              = 1.e-7;
78116d9e3a6SLisandro Dalcin     jac->truncfactor      = 0.0;
78216d9e3a6SLisandro Dalcin     jac->strongthreshold  = .25;
78316d9e3a6SLisandro Dalcin     jac->maxrowsum        = .9;
78416d9e3a6SLisandro Dalcin     jac->coarsentype      = 6;
78516d9e3a6SLisandro Dalcin     jac->measuretype      = 0;
78616d9e3a6SLisandro Dalcin     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = jac->gridsweeps[3]  = 1;
78716d9e3a6SLisandro Dalcin     jac->relaxtype[0]     = jac->relaxtype[1]  = jac->relaxtype[2]  = 3;
78816d9e3a6SLisandro Dalcin     jac->relaxtype[3]     = 9;
78916d9e3a6SLisandro Dalcin     jac->relaxweight      = 1.0;
79016d9e3a6SLisandro Dalcin     jac->outerrelaxweight = 1.0;
79116d9e3a6SLisandro Dalcin     jac->relaxorder       = 1;
79216d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr);
79316d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr);
79416d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr);
79516d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr);
79616d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr);
79716d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr);
79816d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr);
79916d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr);
80016d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr);
80116d9e3a6SLisandro Dalcin     ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr);
80216d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
80316d9e3a6SLisandro Dalcin   }
80416d9e3a6SLisandro Dalcin   ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr);
80516d9e3a6SLisandro Dalcin   jac->hypre_type = PETSC_NULL;
80616d9e3a6SLisandro Dalcin   SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
80716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
80816d9e3a6SLisandro Dalcin }
80916d9e3a6SLisandro Dalcin EXTERN_C_END
81016d9e3a6SLisandro Dalcin 
81116d9e3a6SLisandro Dalcin /*
81216d9e3a6SLisandro Dalcin     It only gets here if the HYPRE type has not been set before the call to
81316d9e3a6SLisandro Dalcin    ...SetFromOptions() which actually is most of the time
81416d9e3a6SLisandro Dalcin */
81516d9e3a6SLisandro Dalcin #undef __FUNCT__
81616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE"
81716d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE(PC pc)
81816d9e3a6SLisandro Dalcin {
81916d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
82016d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
82116d9e3a6SLisandro Dalcin   int            indx;
82216d9e3a6SLisandro Dalcin   const char     *type[] = {"pilut","parasails","boomeramg","euclid"};
82316d9e3a6SLisandro Dalcin   PetscTruth     flg;
82416d9e3a6SLisandro Dalcin 
82516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
82616d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr);
82716d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"pilut",&indx,&flg);CHKERRQ(ierr);
82816d9e3a6SLisandro Dalcin   if (PetscOptionsPublishCount) {   /* force the default if it was not yet set and user did not set with option */
82916d9e3a6SLisandro Dalcin     if (!flg && !jac->hypre_type) {
83016d9e3a6SLisandro Dalcin       flg   = PETSC_TRUE;
83116d9e3a6SLisandro Dalcin       indx = 0;
83216d9e3a6SLisandro Dalcin     }
83316d9e3a6SLisandro Dalcin   }
83416d9e3a6SLisandro Dalcin   if (flg) {
83516d9e3a6SLisandro Dalcin     ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr);
83616d9e3a6SLisandro Dalcin   }
83716d9e3a6SLisandro Dalcin   if (pc->ops->setfromoptions) {
83816d9e3a6SLisandro Dalcin     ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr);
83916d9e3a6SLisandro Dalcin   }
84016d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
84116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
84216d9e3a6SLisandro Dalcin }
84316d9e3a6SLisandro Dalcin 
84416d9e3a6SLisandro Dalcin #undef __FUNCT__
84516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType"
84616d9e3a6SLisandro Dalcin /*@C
84716d9e3a6SLisandro Dalcin      PCHYPRESetType - Sets which hypre preconditioner you wish to use
84816d9e3a6SLisandro Dalcin 
84916d9e3a6SLisandro Dalcin    Input Parameters:
85016d9e3a6SLisandro Dalcin +     pc - the preconditioner context
85116d9e3a6SLisandro Dalcin -     name - either  pilut, parasails, boomeramg, euclid
85216d9e3a6SLisandro Dalcin 
85316d9e3a6SLisandro Dalcin    Options Database Keys:
85416d9e3a6SLisandro Dalcin    -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
85516d9e3a6SLisandro Dalcin 
85616d9e3a6SLisandro Dalcin    Level: intermediate
85716d9e3a6SLisandro Dalcin 
85816d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
85916d9e3a6SLisandro Dalcin            PCHYPRE
86016d9e3a6SLisandro Dalcin 
86116d9e3a6SLisandro Dalcin @*/
86216d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC pc,const char name[])
86316d9e3a6SLisandro Dalcin {
86416d9e3a6SLisandro Dalcin   PetscErrorCode ierr,(*f)(PC,const char[]);
86516d9e3a6SLisandro Dalcin 
86616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
86716d9e3a6SLisandro Dalcin   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
86816d9e3a6SLisandro Dalcin   PetscValidCharPointer(name,2);
86916d9e3a6SLisandro Dalcin   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPRESetType_C",(void (**)(void))&f);CHKERRQ(ierr);
87016d9e3a6SLisandro Dalcin   if (f) {
87116d9e3a6SLisandro Dalcin     ierr = (*f)(pc,name);CHKERRQ(ierr);
87216d9e3a6SLisandro Dalcin   }
87316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
87416d9e3a6SLisandro Dalcin }
87516d9e3a6SLisandro Dalcin 
87616d9e3a6SLisandro Dalcin #undef __FUNCT__
87716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType"
87816d9e3a6SLisandro Dalcin /*@C
87916d9e3a6SLisandro Dalcin      PCHYPREGetType - Gets which hypre preconditioner you are using
88016d9e3a6SLisandro Dalcin 
88116d9e3a6SLisandro Dalcin    Input Parameter:
88216d9e3a6SLisandro Dalcin .     pc - the preconditioner context
88316d9e3a6SLisandro Dalcin 
88416d9e3a6SLisandro Dalcin    Output Parameter:
88516d9e3a6SLisandro Dalcin .     name - either  pilut, parasails, boomeramg, euclid
88616d9e3a6SLisandro Dalcin 
88716d9e3a6SLisandro Dalcin    Level: intermediate
88816d9e3a6SLisandro Dalcin 
88916d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC,
89016d9e3a6SLisandro Dalcin            PCHYPRE
89116d9e3a6SLisandro Dalcin 
89216d9e3a6SLisandro Dalcin @*/
89316d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC pc,const char *name[])
89416d9e3a6SLisandro Dalcin {
89516d9e3a6SLisandro Dalcin   PetscErrorCode ierr,(*f)(PC,const char*[]);
89616d9e3a6SLisandro Dalcin 
89716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
89816d9e3a6SLisandro Dalcin   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
89916d9e3a6SLisandro Dalcin   PetscValidPointer(name,2);
90016d9e3a6SLisandro Dalcin   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPREGetType_C",(void (**)(void))&f);CHKERRQ(ierr);
90116d9e3a6SLisandro Dalcin   if (f) {
90216d9e3a6SLisandro Dalcin     ierr = (*f)(pc,name);CHKERRQ(ierr);
90316d9e3a6SLisandro Dalcin   }
90416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
90516d9e3a6SLisandro Dalcin }
90616d9e3a6SLisandro Dalcin 
90716d9e3a6SLisandro Dalcin /*MC
90816d9e3a6SLisandro Dalcin      PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre
90916d9e3a6SLisandro Dalcin 
91016d9e3a6SLisandro Dalcin    Options Database Keys:
91116d9e3a6SLisandro Dalcin +   -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
91216d9e3a6SLisandro Dalcin -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
91316d9e3a6SLisandro Dalcin           preconditioner
91416d9e3a6SLisandro Dalcin 
91516d9e3a6SLisandro Dalcin    Level: intermediate
91616d9e3a6SLisandro Dalcin 
91716d9e3a6SLisandro Dalcin    Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()),
91816d9e3a6SLisandro Dalcin           the many hypre options can ONLY be set via the options database (e.g. the command line
91916d9e3a6SLisandro Dalcin           or with PetscOptionsSetValue(), there are no functions to set them)
92016d9e3a6SLisandro Dalcin 
92116d9e3a6SLisandro Dalcin           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations
92216d9e3a6SLisandro Dalcin           (V-cycles) that boomeramg does EACH time it is called. So for example, if it is set to 2 then
92316d9e3a6SLisandro Dalcin           2-V-cycles are being used to define the preconditioner. -ksp_max_iter and -ksp_rtol STILL determine
92416d9e3a6SLisandro Dalcin           the total number of iterations and tolerance for the Krylov solver. For example, if
92516d9e3a6SLisandro Dalcin           -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10 then AT MOST twenty V-cycles of boomeramg
92616d9e3a6SLisandro Dalcin           will be called.
92716d9e3a6SLisandro Dalcin 
92816d9e3a6SLisandro Dalcin           If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly
92916d9e3a6SLisandro Dalcin           and use -ksp_max_it to control the number of V-cycles.
93016d9e3a6SLisandro Dalcin           (see the PETSc FAQ.html at the PETSc website under the Documentation tab).
93116d9e3a6SLisandro Dalcin 
93216d9e3a6SLisandro Dalcin 	  2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option
93316d9e3a6SLisandro Dalcin 	  -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L.
93416d9e3a6SLisandro Dalcin 
93516d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
93616d9e3a6SLisandro Dalcin            PCHYPRESetType()
93716d9e3a6SLisandro Dalcin 
93816d9e3a6SLisandro Dalcin M*/
93916d9e3a6SLisandro Dalcin 
94016d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
94116d9e3a6SLisandro Dalcin #undef __FUNCT__
94216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCCreate_HYPRE"
94316d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_HYPRE(PC pc)
94416d9e3a6SLisandro Dalcin {
94516d9e3a6SLisandro Dalcin   PC_HYPRE       *jac;
94616d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
94716d9e3a6SLisandro Dalcin 
94816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
94916d9e3a6SLisandro Dalcin   ierr = PetscNew(PC_HYPRE,&jac);CHKERRQ(ierr);
95016d9e3a6SLisandro Dalcin   ierr = PetscLogObjectMemory(pc,sizeof(PC_HYPRE));CHKERRQ(ierr);
95116d9e3a6SLisandro Dalcin   pc->data                 = jac;
95216d9e3a6SLisandro Dalcin   pc->ops->destroy         = PCDestroy_HYPRE;
95316d9e3a6SLisandro Dalcin   pc->ops->setfromoptions  = PCSetFromOptions_HYPRE;
95416d9e3a6SLisandro Dalcin   pc->ops->setup           = PCSetUp_HYPRE;
95516d9e3a6SLisandro Dalcin   pc->ops->apply           = PCApply_HYPRE;
95616d9e3a6SLisandro Dalcin   jac->comm_hypre          = MPI_COMM_NULL;
95716d9e3a6SLisandro Dalcin   jac->hypre_type          = PETSC_NULL;
95816d9e3a6SLisandro Dalcin   /* duplicate communicator for hypre */
95916d9e3a6SLisandro Dalcin   ierr = MPI_Comm_dup(pc->comm,&(jac->comm_hypre));CHKERRQ(ierr);
96016d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr);
96116d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr);
96216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
96316d9e3a6SLisandro Dalcin }
96416d9e3a6SLisandro Dalcin EXTERN_C_END
965