xref: /petsc/src/ksp/pc/impls/hypre/hypre.c (revision b00a91154f763f12aa55f3d53a3f2776f15f49e3)
116d9e3a6SLisandro Dalcin 
216d9e3a6SLisandro Dalcin /*
316d9e3a6SLisandro Dalcin    Provides an interface to the LLNL package hypre
416d9e3a6SLisandro Dalcin */
50f1074feSSatish Balay 
60f1074feSSatish Balay /* Must use hypre 2.0.0 or more recent. */
70f1074feSSatish Balay 
8b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>          /*I "petscpc.h" I*/
9c6db04a5SJed Brown #include <../src/dm/impls/da/hypre/mhyp.h>
1016d9e3a6SLisandro Dalcin 
11dff31646SBarry Smith static PetscBool cite = PETSC_FALSE;
121f817a21SBarry Smith static const char hypreCitation[] = "@manual{hypre-web-page,\n  title  = {{\\sl hypre}: High Performance Preconditioners},\n  organization = {Lawrence Livermore National Laboratory},\n  note  = {\\url{http://www.llnl.gov/CASC/hypre/}}\n}\n";
131f817a21SBarry Smith 
1416d9e3a6SLisandro Dalcin /*
1516d9e3a6SLisandro Dalcin    Private context (data structure) for the  preconditioner.
1616d9e3a6SLisandro Dalcin */
1716d9e3a6SLisandro Dalcin typedef struct {
1816d9e3a6SLisandro Dalcin   HYPRE_Solver   hsolver;
1916d9e3a6SLisandro Dalcin   HYPRE_IJMatrix ij;
2016d9e3a6SLisandro Dalcin   HYPRE_IJVector b,x;
2116d9e3a6SLisandro Dalcin 
224ddd07fcSJed Brown   HYPRE_Int (*destroy)(HYPRE_Solver);
234ddd07fcSJed Brown   HYPRE_Int (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
244ddd07fcSJed Brown   HYPRE_Int (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
2516d9e3a6SLisandro Dalcin 
2616d9e3a6SLisandro Dalcin   MPI_Comm comm_hypre;
2716d9e3a6SLisandro Dalcin   char     *hypre_type;
2816d9e3a6SLisandro Dalcin 
2916d9e3a6SLisandro Dalcin   /* options for Pilut and BoomerAMG*/
304ddd07fcSJed Brown   PetscInt maxiter;
3116d9e3a6SLisandro Dalcin   double   tol;
3216d9e3a6SLisandro Dalcin 
3316d9e3a6SLisandro Dalcin   /* options for Pilut */
344ddd07fcSJed Brown   PetscInt factorrowsize;
3516d9e3a6SLisandro Dalcin 
3616d9e3a6SLisandro Dalcin   /* options for ParaSails */
374ddd07fcSJed Brown   PetscInt nlevels;
3816d9e3a6SLisandro Dalcin   double   threshhold;
3916d9e3a6SLisandro Dalcin   double   filter;
404ddd07fcSJed Brown   PetscInt sym;
4116d9e3a6SLisandro Dalcin   double   loadbal;
424ddd07fcSJed Brown   PetscInt logging;
434ddd07fcSJed Brown   PetscInt ruse;
444ddd07fcSJed Brown   PetscInt symt;
4516d9e3a6SLisandro Dalcin 
4616d9e3a6SLisandro Dalcin   /* options for Euclid */
47ace3abfcSBarry Smith   PetscBool bjilu;
484ddd07fcSJed Brown   PetscInt  levels;
4916d9e3a6SLisandro Dalcin 
5016d9e3a6SLisandro Dalcin   /* options for Euclid and BoomerAMG */
51ace3abfcSBarry Smith   PetscBool printstatistics;
5216d9e3a6SLisandro Dalcin 
5316d9e3a6SLisandro Dalcin   /* options for BoomerAMG */
544ddd07fcSJed Brown   PetscInt  cycletype;
554ddd07fcSJed Brown   PetscInt  maxlevels;
5616d9e3a6SLisandro Dalcin   double    strongthreshold;
5716d9e3a6SLisandro Dalcin   double    maxrowsum;
584ddd07fcSJed Brown   PetscInt  gridsweeps[3];
594ddd07fcSJed Brown   PetscInt  coarsentype;
604ddd07fcSJed Brown   PetscInt  measuretype;
614ddd07fcSJed Brown   PetscInt  relaxtype[3];
6216d9e3a6SLisandro Dalcin   double    relaxweight;
6316d9e3a6SLisandro Dalcin   double    outerrelaxweight;
644ddd07fcSJed Brown   PetscInt  relaxorder;
6516d9e3a6SLisandro Dalcin   double    truncfactor;
66ace3abfcSBarry Smith   PetscBool applyrichardson;
674ddd07fcSJed Brown   PetscInt  pmax;
684ddd07fcSJed Brown   PetscInt  interptype;
694ddd07fcSJed Brown   PetscInt  agg_nl;
704ddd07fcSJed Brown   PetscInt  agg_num_paths;
714ddd07fcSJed Brown   PetscInt  nodal_coarsen;
72ace3abfcSBarry Smith   PetscBool nodal_relax;
734ddd07fcSJed Brown   PetscInt  nodal_relax_levels;
7416d9e3a6SLisandro Dalcin } PC_HYPRE;
7516d9e3a6SLisandro Dalcin 
76d2128fa2SBarry Smith #undef __FUNCT__
77d2128fa2SBarry Smith #define __FUNCT__ "PCHYPREGetSolver"
78d2128fa2SBarry Smith PetscErrorCode PCHYPREGetSolver(PC pc,HYPRE_Solver *hsolver)
79d2128fa2SBarry Smith {
80d2128fa2SBarry Smith   PC_HYPRE *jac = (PC_HYPRE*)pc->data;
81d2128fa2SBarry Smith 
82d2128fa2SBarry Smith   PetscFunctionBegin;
83d2128fa2SBarry Smith   *hsolver = jac->hsolver;
84d2128fa2SBarry Smith   PetscFunctionReturn(0);
85d2128fa2SBarry Smith }
8616d9e3a6SLisandro Dalcin 
8716d9e3a6SLisandro Dalcin #undef __FUNCT__
8816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetUp_HYPRE"
8916d9e3a6SLisandro Dalcin static PetscErrorCode PCSetUp_HYPRE(PC pc)
9016d9e3a6SLisandro Dalcin {
9116d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
9216d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
9316d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
9416d9e3a6SLisandro Dalcin   HYPRE_ParVector    bv,xv;
9516d9e3a6SLisandro Dalcin   PetscInt           bs;
9616d9e3a6SLisandro Dalcin 
9716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9816d9e3a6SLisandro Dalcin   if (!jac->hypre_type) {
9902a17cd4SBarry Smith     ierr = PCHYPRESetType(pc,"boomeramg");CHKERRQ(ierr);
10016d9e3a6SLisandro Dalcin   }
1015f5c5b43SBarry Smith 
1025f5c5b43SBarry Smith   if (pc->setupcalled) {
1035f5c5b43SBarry Smith     /* always destroy the old matrix and create a new memory;
1045f5c5b43SBarry Smith        hope this does not churn the memory too much. The problem
1055f5c5b43SBarry Smith        is I do not know if it is possible to put the matrix back to
1065f5c5b43SBarry Smith        its initial state so that we can directly copy the values
1075f5c5b43SBarry Smith        the second time through. */
108fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->ij));
1095f5c5b43SBarry Smith     jac->ij = 0;
11016d9e3a6SLisandro Dalcin   }
1115f5c5b43SBarry Smith 
11216d9e3a6SLisandro Dalcin   if (!jac->ij) { /* create the matrix the first time through */
11316d9e3a6SLisandro Dalcin     ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr);
11416d9e3a6SLisandro Dalcin   }
11516d9e3a6SLisandro Dalcin   if (!jac->b) { /* create the vectors the first time through */
11616d9e3a6SLisandro Dalcin     Vec x,b;
11716d9e3a6SLisandro Dalcin     ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr);
11816d9e3a6SLisandro Dalcin     ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr);
11916d9e3a6SLisandro Dalcin     ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr);
1206bf464f9SBarry Smith     ierr = VecDestroy(&x);CHKERRQ(ierr);
1216bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
12216d9e3a6SLisandro Dalcin   }
1235f5c5b43SBarry Smith 
12416d9e3a6SLisandro Dalcin   /* special case for BoomerAMG */
12516d9e3a6SLisandro Dalcin   if (jac->setup == HYPRE_BoomerAMGSetup) {
12616d9e3a6SLisandro Dalcin     ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr);
1272fa5cd67SKarl Rupp     if (bs > 1) PetscStackCallStandard(HYPRE_BoomerAMGSetNumFunctions,(jac->hsolver,bs));
12816d9e3a6SLisandro Dalcin   };
12916d9e3a6SLisandro Dalcin   ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr);
130fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
131fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&bv));
132fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&xv));
133fd3f9acdSBarry Smith   PetscStackCall("HYPRE_SetupXXX",ierr = (*jac->setup)(jac->hsolver,hmat,bv,xv);CHKERRQ(ierr););
13416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
13516d9e3a6SLisandro Dalcin }
13616d9e3a6SLisandro Dalcin 
13716d9e3a6SLisandro Dalcin /*
13816d9e3a6SLisandro Dalcin     Replaces the address where the HYPRE vector points to its data with the address of
13916d9e3a6SLisandro Dalcin   PETSc's data. Saves the old address so it can be reset when we are finished with it.
14016d9e3a6SLisandro Dalcin   Allows use to get the data into a HYPRE vector without the cost of memcopies
14116d9e3a6SLisandro Dalcin */
14216d9e3a6SLisandro Dalcin #define HYPREReplacePointer(b,newvalue,savedvalue) { \
14316d9e3a6SLisandro Dalcin     hypre_ParVector *par_vector   = (hypre_ParVector*)hypre_IJVectorObject(((hypre_IJVector*)b)); \
14416d9e3a6SLisandro Dalcin     hypre_Vector    *local_vector = hypre_ParVectorLocalVector(par_vector); \
14516d9e3a6SLisandro Dalcin     savedvalue         = local_vector->data; \
1460ad7597dSKarl Rupp     local_vector->data = newvalue;          \
1470ad7597dSKarl Rupp }
14816d9e3a6SLisandro Dalcin 
14916d9e3a6SLisandro Dalcin #undef __FUNCT__
15016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApply_HYPRE"
15116d9e3a6SLisandro Dalcin static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
15216d9e3a6SLisandro Dalcin {
15316d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
15416d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
15516d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
15616d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
15716d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
15816d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
1594ddd07fcSJed Brown   PetscInt           hierr;
16016d9e3a6SLisandro Dalcin 
16116d9e3a6SLisandro Dalcin   PetscFunctionBegin;
162dff31646SBarry Smith   ierr = PetscCitationsRegister(hypreCitation,&cite);CHKERRQ(ierr);
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 
169fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
170fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&jbv));
171fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&jxv));
172fd3f9acdSBarry Smith   PetscStackCall("Hypre solve",hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv);
17365e19b50SBarry Smith                                if (hierr && hierr != HYPRE_ERROR_CONV) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
174fd3f9acdSBarry Smith                                if (hierr) hypre__global_error = 0;);
17516d9e3a6SLisandro Dalcin 
17616d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
17716d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
17816d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
17916d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
18016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
18116d9e3a6SLisandro Dalcin }
18216d9e3a6SLisandro Dalcin 
18316d9e3a6SLisandro Dalcin #undef __FUNCT__
18416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCDestroy_HYPRE"
18516d9e3a6SLisandro Dalcin static PetscErrorCode PCDestroy_HYPRE(PC pc)
18616d9e3a6SLisandro Dalcin {
18716d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
18816d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
18916d9e3a6SLisandro Dalcin 
19016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
191fd3f9acdSBarry Smith   if (jac->ij) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->ij));
192fd3f9acdSBarry Smith   if (jac->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->b));
193fd3f9acdSBarry Smith   if (jac->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->x));
194226b0620SJed Brown   if (jac->destroy) PetscStackCall("HYPRE_DestroyXXX",ierr = (*jac->destroy)(jac->hsolver);CHKERRQ(ierr););
195503cfb0cSBarry Smith   ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr);
19616d9e3a6SLisandro Dalcin   if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);}
197c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
19816d9e3a6SLisandro Dalcin 
19916d9e3a6SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr);
200bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetType_C",NULL);CHKERRQ(ierr);
201bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHYPREGetType_C",NULL);CHKERRQ(ierr);
20216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
20316d9e3a6SLisandro Dalcin }
20416d9e3a6SLisandro Dalcin 
20516d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
20616d9e3a6SLisandro Dalcin #undef __FUNCT__
20716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut"
20816d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc)
20916d9e3a6SLisandro Dalcin {
21016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
21116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
212ace3abfcSBarry Smith   PetscBool      flag;
21316d9e3a6SLisandro Dalcin 
21416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
21516d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr);
21616d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr);
217fd3f9acdSBarry Smith   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetMaxIter,(jac->hsolver,jac->maxiter));
21816d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr);
219fd3f9acdSBarry Smith   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetDropTolerance,(jac->hsolver,jac->tol));
22016d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr);
221fd3f9acdSBarry Smith   if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetFactorRowSize,(jac->hsolver,jac->factorrowsize));
22216d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
22316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
22416d9e3a6SLisandro Dalcin }
22516d9e3a6SLisandro Dalcin 
22616d9e3a6SLisandro Dalcin #undef __FUNCT__
22716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Pilut"
22816d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
22916d9e3a6SLisandro Dalcin {
23016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
23116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
232ace3abfcSBarry Smith   PetscBool      iascii;
23316d9e3a6SLisandro Dalcin 
23416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
235251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
23616d9e3a6SLisandro Dalcin   if (iascii) {
23716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");CHKERRQ(ierr);
23816d9e3a6SLisandro Dalcin     if (jac->maxiter != PETSC_DEFAULT) {
23916d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr);
24016d9e3a6SLisandro Dalcin     } else {
24116d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr);
24216d9e3a6SLisandro Dalcin     }
24316d9e3a6SLisandro Dalcin     if (jac->tol != PETSC_DEFAULT) {
24416d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr);
24516d9e3a6SLisandro Dalcin     } else {
24616d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr);
24716d9e3a6SLisandro Dalcin     }
24816d9e3a6SLisandro Dalcin     if (jac->factorrowsize != PETSC_DEFAULT) {
24916d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr);
25016d9e3a6SLisandro Dalcin     } else {
25116d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");CHKERRQ(ierr);
25216d9e3a6SLisandro Dalcin     }
25316d9e3a6SLisandro Dalcin   }
25416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
25516d9e3a6SLisandro Dalcin }
25616d9e3a6SLisandro Dalcin 
25716d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
25816d9e3a6SLisandro Dalcin #undef __FUNCT__
25916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid"
26016d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc)
26116d9e3a6SLisandro Dalcin {
26216d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
26316d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
264ace3abfcSBarry Smith   PetscBool      flag;
265390e7148SBarry Smith   char           *args[8],levels[16];
266390e7148SBarry Smith   PetscInt       cnt = 0;
26716d9e3a6SLisandro Dalcin 
26816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
26916d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr);
27016d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr);
27116d9e3a6SLisandro Dalcin   if (flag) {
272ce94432eSBarry Smith     if (jac->levels < 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels);
2738caf3d72SBarry Smith     ierr        = PetscSNPrintf(levels,sizeof(levels),"%D",jac->levels);CHKERRQ(ierr);
274390e7148SBarry Smith     args[cnt++] = (char*)"-level"; args[cnt++] = levels;
27516d9e3a6SLisandro Dalcin   }
2760298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,NULL);CHKERRQ(ierr);
27716d9e3a6SLisandro Dalcin   if (jac->bjilu) {
278390e7148SBarry Smith     args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1";
27916d9e3a6SLisandro Dalcin   }
28016d9e3a6SLisandro Dalcin 
2810298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,NULL);CHKERRQ(ierr);
28216d9e3a6SLisandro Dalcin   if (jac->printstatistics) {
283390e7148SBarry Smith     args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1";
284390e7148SBarry Smith     args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1";
28516d9e3a6SLisandro Dalcin   }
28616d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
287fd3f9acdSBarry Smith   if (cnt) PetscStackCallStandard(HYPRE_EuclidSetParams,(jac->hsolver,cnt,args));
28816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
28916d9e3a6SLisandro Dalcin }
29016d9e3a6SLisandro Dalcin 
29116d9e3a6SLisandro Dalcin #undef __FUNCT__
29216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Euclid"
29316d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
29416d9e3a6SLisandro Dalcin {
29516d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
29616d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
297ace3abfcSBarry Smith   PetscBool      iascii;
29816d9e3a6SLisandro Dalcin 
29916d9e3a6SLisandro Dalcin   PetscFunctionBegin;
300251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
30116d9e3a6SLisandro Dalcin   if (iascii) {
30216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n");CHKERRQ(ierr);
30316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr);
30416d9e3a6SLisandro Dalcin     if (jac->bjilu) {
30516d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr);
30616d9e3a6SLisandro Dalcin     }
30716d9e3a6SLisandro Dalcin   }
30816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
30916d9e3a6SLisandro Dalcin }
31016d9e3a6SLisandro Dalcin 
31116d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
31216d9e3a6SLisandro Dalcin 
31316d9e3a6SLisandro Dalcin #undef __FUNCT__
31416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG"
31516d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
31616d9e3a6SLisandro Dalcin {
31716d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
31816d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
31916d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
32016d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
32116d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
32216d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
3234ddd07fcSJed Brown   PetscInt           hierr;
32416d9e3a6SLisandro Dalcin 
32516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
326dff31646SBarry Smith   ierr = PetscCitationsRegister(hypreCitation,&cite);CHKERRQ(ierr);
32716d9e3a6SLisandro Dalcin   ierr = VecSet(x,0.0);CHKERRQ(ierr);
32816d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
32916d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
33016d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
33116d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
33216d9e3a6SLisandro Dalcin 
333fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
334fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&jbv));
335fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&jxv));
33616d9e3a6SLisandro Dalcin 
33716d9e3a6SLisandro Dalcin   hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv);
33816d9e3a6SLisandro Dalcin   /* error code of 1 in BoomerAMG merely means convergence not achieved */
339e32f2f54SBarry Smith   if (hierr && (hierr != 1)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
34016d9e3a6SLisandro Dalcin   if (hierr) hypre__global_error = 0;
34116d9e3a6SLisandro Dalcin 
34216d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
34316d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
34416d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
34516d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
34616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
34716d9e3a6SLisandro Dalcin }
34816d9e3a6SLisandro Dalcin 
349a669f990SJed Brown /* static array length */
350a669f990SJed Brown #define ALEN(a) (sizeof(a)/sizeof((a)[0]))
351a669f990SJed Brown 
35216d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
3530f1074feSSatish Balay static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"};
35416d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
35565de4495SJed Brown /* The following corresponds to HYPRE_BoomerAMGSetRelaxType which has many missing numbers in the enum */
35665de4495SJed Brown static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","seqboundary-Gauss-Seidel","SOR/Jacobi","backward-SOR/Jacobi",
35765de4495SJed Brown                                                   "" /* [5] hybrid chaotic Gauss-Seidel (works only with OpenMP) */,"symmetric-SOR/Jacobi",
35865de4495SJed Brown                                                   "" /* 7 */,"l1scaled-SOR/Jacobi","Gaussian-elimination",
35965de4495SJed Brown                                                   "" /* 10 */, "" /* 11 */, "" /* 12 */, "" /* 13 */, "" /* 14 */,
36065de4495SJed Brown                                                   "CG" /* non-stationary */,"Chebyshev","FCF-Jacobi","l1scaled-Jacobi"};
3610f1074feSSatish Balay static const char *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i",
3620f1074feSSatish Balay                                                   "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"};
36316d9e3a6SLisandro Dalcin #undef __FUNCT__
36416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG"
36516d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
36616d9e3a6SLisandro Dalcin {
36716d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
36816d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
3694ddd07fcSJed Brown   PetscInt       n,indx,level;
370ace3abfcSBarry Smith   PetscBool      flg, tmp_truth;
37116d9e3a6SLisandro Dalcin   double         tmpdbl, twodbl[2];
37216d9e3a6SLisandro Dalcin 
37316d9e3a6SLisandro Dalcin   PetscFunctionBegin;
37416d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr);
3754336a9eeSBarry Smith   ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType+1,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr);
37616d9e3a6SLisandro Dalcin   if (flg) {
3774336a9eeSBarry Smith     jac->cycletype = indx+1;
378fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
37916d9e3a6SLisandro Dalcin   }
38016d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr);
38116d9e3a6SLisandro Dalcin   if (flg) {
382ce94432eSBarry Smith     if (jac->maxlevels < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
383fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
38416d9e3a6SLisandro Dalcin   }
38516d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr);
38616d9e3a6SLisandro Dalcin   if (flg) {
387ce94432eSBarry Smith     if (jac->maxiter < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
388fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
38916d9e3a6SLisandro Dalcin   }
3900f1074feSSatish 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);
39116d9e3a6SLisandro Dalcin   if (flg) {
392ce94432eSBarry Smith     if (jac->tol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be greater than or equal to zero",jac->tol);
393fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
39416d9e3a6SLisandro Dalcin   }
39516d9e3a6SLisandro Dalcin 
3960f1074feSSatish Balay   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr);
39716d9e3a6SLisandro Dalcin   if (flg) {
398ce94432eSBarry Smith     if (jac->truncfactor < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor);
399fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
40016d9e3a6SLisandro Dalcin   }
40116d9e3a6SLisandro Dalcin 
4020f1074feSSatish 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);
4030f1074feSSatish Balay   if (flg) {
404ce94432eSBarry Smith     if (jac->pmax < 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"P_max %G must be greater than or equal to zero",jac->pmax);
405fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
4060f1074feSSatish Balay   }
4070f1074feSSatish Balay 
4080f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);CHKERRQ(ierr);
4090f1074feSSatish Balay   if (flg) {
410ce94432eSBarry Smith     if (jac->agg_nl < 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %G must be greater than or equal to zero",jac->agg_nl);
4110f1074feSSatish Balay 
412fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
4130f1074feSSatish Balay   }
4140f1074feSSatish Balay 
4150f1074feSSatish Balay 
4160f1074feSSatish 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);
4170f1074feSSatish Balay   if (flg) {
418ce94432eSBarry Smith     if (jac->agg_num_paths < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of paths %G must be greater than or equal to 1",jac->agg_num_paths);
4190f1074feSSatish Balay 
420fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
4210f1074feSSatish Balay   }
4220f1074feSSatish Balay 
4230f1074feSSatish Balay 
42416d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr);
42516d9e3a6SLisandro Dalcin   if (flg) {
426ce94432eSBarry Smith     if (jac->strongthreshold < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold);
427fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
42816d9e3a6SLisandro Dalcin   }
42916d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr);
43016d9e3a6SLisandro Dalcin   if (flg) {
431ce94432eSBarry Smith     if (jac->maxrowsum < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum);
432ce94432eSBarry Smith     if (jac->maxrowsum > 1.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum);
433fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
43416d9e3a6SLisandro Dalcin   }
43516d9e3a6SLisandro Dalcin 
43616d9e3a6SLisandro Dalcin   /* Grid sweeps */
4370f1074feSSatish 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);
43816d9e3a6SLisandro Dalcin   if (flg) {
439fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver,indx));
44016d9e3a6SLisandro Dalcin     /* modify the jac structure so we can view the updated options with PC_View */
44116d9e3a6SLisandro Dalcin     jac->gridsweeps[0] = indx;
4420f1074feSSatish Balay     jac->gridsweeps[1] = indx;
4430f1074feSSatish Balay     /*defaults coarse to 1 */
4440f1074feSSatish Balay     jac->gridsweeps[2] = 1;
44516d9e3a6SLisandro Dalcin   }
4460f1074feSSatish Balay 
4470f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx,&flg);CHKERRQ(ierr);
44816d9e3a6SLisandro Dalcin   if (flg) {
449fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 1));
4500f1074feSSatish Balay     jac->gridsweeps[0] = indx;
45116d9e3a6SLisandro Dalcin   }
45216d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr);
45316d9e3a6SLisandro Dalcin   if (flg) {
454fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 2));
4550f1074feSSatish Balay     jac->gridsweeps[1] = indx;
45616d9e3a6SLisandro Dalcin   }
4570f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr);
45816d9e3a6SLisandro Dalcin   if (flg) {
459fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 3));
4600f1074feSSatish Balay     jac->gridsweeps[2] = indx;
46116d9e3a6SLisandro Dalcin   }
46216d9e3a6SLisandro Dalcin 
46316d9e3a6SLisandro Dalcin   /* Relax type */
464a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_all","Relax type for the up and down cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
46516d9e3a6SLisandro Dalcin   if (flg) {
4660f1074feSSatish Balay     jac->relaxtype[0] = jac->relaxtype[1]  = indx;
467fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, indx));
4680f1074feSSatish Balay     /* by default, coarse type set to 9 */
4690f1074feSSatish Balay     jac->relaxtype[2] = 9;
4700f1074feSSatish Balay 
47116d9e3a6SLisandro Dalcin   }
472a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
47316d9e3a6SLisandro Dalcin   if (flg) {
47416d9e3a6SLisandro Dalcin     jac->relaxtype[0] = indx;
475fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 1));
47616d9e3a6SLisandro Dalcin   }
477a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
47816d9e3a6SLisandro Dalcin   if (flg) {
4790f1074feSSatish Balay     jac->relaxtype[1] = indx;
480fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 2));
48116d9e3a6SLisandro Dalcin   }
482a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr);
48316d9e3a6SLisandro Dalcin   if (flg) {
4840f1074feSSatish Balay     jac->relaxtype[2] = indx;
485fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 3));
48616d9e3a6SLisandro Dalcin   }
48716d9e3a6SLisandro Dalcin 
48816d9e3a6SLisandro Dalcin   /* Relaxation Weight */
48916d9e3a6SLisandro 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);
49016d9e3a6SLisandro Dalcin   if (flg) {
491fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxWt,(jac->hsolver,tmpdbl));
49216d9e3a6SLisandro Dalcin     jac->relaxweight = tmpdbl;
49316d9e3a6SLisandro Dalcin   }
49416d9e3a6SLisandro Dalcin 
49516d9e3a6SLisandro Dalcin   n         = 2;
49616d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
49716d9e3a6SLisandro 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);
49816d9e3a6SLisandro Dalcin   if (flg) {
49916d9e3a6SLisandro Dalcin     if (n == 2) {
50016d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
501fd3f9acdSBarry Smith       PetscStackCallStandard(HYPRE_BoomerAMGSetLevelRelaxWt,(jac->hsolver,twodbl[0],indx));
502ce94432eSBarry Smith     } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relax weight level: you must provide 2 values separated by a comma (and no space), you provided %d",n);
50316d9e3a6SLisandro Dalcin   }
50416d9e3a6SLisandro Dalcin 
50516d9e3a6SLisandro Dalcin   /* Outer relaxation Weight */
50616d9e3a6SLisandro 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);
50716d9e3a6SLisandro Dalcin   if (flg) {
508fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetOuterWt,(jac->hsolver, tmpdbl));
50916d9e3a6SLisandro Dalcin     jac->outerrelaxweight = tmpdbl;
51016d9e3a6SLisandro Dalcin   }
51116d9e3a6SLisandro Dalcin 
51216d9e3a6SLisandro Dalcin   n         = 2;
51316d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
51416d9e3a6SLisandro 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);
51516d9e3a6SLisandro Dalcin   if (flg) {
51616d9e3a6SLisandro Dalcin     if (n == 2) {
51716d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
518fd3f9acdSBarry Smith       PetscStackCallStandard(HYPRE_BoomerAMGSetLevelOuterWt,(jac->hsolver, twodbl[0], indx));
519ce94432eSBarry Smith     } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relax weight outer level: You must provide 2 values separated by a comma (and no space), you provided %d",n);
52016d9e3a6SLisandro Dalcin   }
52116d9e3a6SLisandro Dalcin 
52216d9e3a6SLisandro Dalcin   /* the Relax Order */
523acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
52416d9e3a6SLisandro Dalcin 
52516d9e3a6SLisandro Dalcin   if (flg) {
52616d9e3a6SLisandro Dalcin     jac->relaxorder = 0;
527fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
52816d9e3a6SLisandro Dalcin   }
529a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,ALEN(HYPREBoomerAMGMeasureType),HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr);
53016d9e3a6SLisandro Dalcin   if (flg) {
53116d9e3a6SLisandro Dalcin     jac->measuretype = indx;
532fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
53316d9e3a6SLisandro Dalcin   }
5340f1074feSSatish Balay   /* update list length 3/07 */
535a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,ALEN(HYPREBoomerAMGCoarsenType),HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr);
53616d9e3a6SLisandro Dalcin   if (flg) {
53716d9e3a6SLisandro Dalcin     jac->coarsentype = indx;
538fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
53916d9e3a6SLisandro Dalcin   }
5400f1074feSSatish Balay 
5410f1074feSSatish Balay   /* new 3/07 */
542a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,ALEN(HYPREBoomerAMGInterpType),HYPREBoomerAMGInterpType[0],&indx,&flg);CHKERRQ(ierr);
5430f1074feSSatish Balay   if (flg) {
5440f1074feSSatish Balay     jac->interptype = indx;
545fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
5460f1074feSSatish Balay   }
5470f1074feSSatish Balay 
548b96a4a96SBarry Smith   ierr = PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);CHKERRQ(ierr);
54916d9e3a6SLisandro Dalcin   if (flg) {
550b96a4a96SBarry Smith     level = 3;
5510298fd71SBarry Smith     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,NULL);CHKERRQ(ierr);
5522fa5cd67SKarl Rupp 
553b96a4a96SBarry Smith     jac->printstatistics = PETSC_TRUE;
554fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetPrintLevel,(jac->hsolver,level));
5552ae77aedSBarry Smith   }
5562ae77aedSBarry Smith 
557b96a4a96SBarry Smith   ierr = PetscOptionsName("-pc_hypre_boomeramg_print_debug","Print debug information","None",&flg);CHKERRQ(ierr);
5582ae77aedSBarry Smith   if (flg) {
559b96a4a96SBarry Smith     level = 3;
5600298fd71SBarry Smith     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,NULL);CHKERRQ(ierr);
5612fa5cd67SKarl Rupp 
562b96a4a96SBarry Smith     jac->printstatistics = PETSC_TRUE;
563fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetDebugFlag,(jac->hsolver,level));
56416d9e3a6SLisandro Dalcin   }
5658f87f92bSBarry Smith 
566acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5678f87f92bSBarry Smith   if (flg && tmp_truth) {
5688f87f92bSBarry Smith     jac->nodal_coarsen = 1;
569fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetNodal,(jac->hsolver,1));
5708f87f92bSBarry Smith   }
5718f87f92bSBarry Smith 
572acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5738f87f92bSBarry Smith   if (flg && tmp_truth) {
5748f87f92bSBarry Smith     PetscInt tmp_int;
5758f87f92bSBarry Smith     ierr = PetscOptionsInt("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);CHKERRQ(ierr);
5768f87f92bSBarry Smith     if (flg) jac->nodal_relax_levels = tmp_int;
577fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothType,(jac->hsolver,6));
578fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetDomainType,(jac->hsolver,1));
579fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetOverlap,(jac->hsolver,0));
580fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothNumLevels,(jac->hsolver,jac->nodal_relax_levels));
5818f87f92bSBarry Smith   }
5828f87f92bSBarry Smith 
58316d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
58416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
58516d9e3a6SLisandro Dalcin }
58616d9e3a6SLisandro Dalcin 
58716d9e3a6SLisandro Dalcin #undef __FUNCT__
58816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG"
589ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
59016d9e3a6SLisandro Dalcin {
59116d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
59216d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
5934ddd07fcSJed Brown   PetscInt       oits;
59416d9e3a6SLisandro Dalcin 
59516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
596dff31646SBarry Smith   ierr = PetscCitationsRegister(hypreCitation,&cite);CHKERRQ(ierr);
597fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,its*jac->maxiter));
598fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,rtol));
59916d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_TRUE;
60016d9e3a6SLisandro Dalcin   ierr                 = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr);
60116d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_FALSE;
6028b1f7689SBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGGetNumIterations,(jac->hsolver,(HYPRE_Int *)&oits));
6034d0a8057SBarry Smith   *outits = oits;
6044d0a8057SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
6054d0a8057SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
606fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
607fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
60816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
60916d9e3a6SLisandro Dalcin }
61016d9e3a6SLisandro Dalcin 
61116d9e3a6SLisandro Dalcin 
61216d9e3a6SLisandro Dalcin #undef __FUNCT__
61316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_BoomerAMG"
61416d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
61516d9e3a6SLisandro Dalcin {
61616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
61716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
618ace3abfcSBarry Smith   PetscBool      iascii;
61916d9e3a6SLisandro Dalcin 
62016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
621251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
62216d9e3a6SLisandro Dalcin   if (iascii) {
62316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr);
62416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr);
62516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr);
62616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr);
62716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr);
62816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr);
6290f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation truncation factor %G\n",jac->truncfactor);CHKERRQ(ierr);
6300f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);CHKERRQ(ierr);
6310f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);CHKERRQ(ierr);
6320f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);CHKERRQ(ierr);
6330f1074feSSatish Balay 
63416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr);
63516d9e3a6SLisandro Dalcin 
6360f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps down         %d\n",jac->gridsweeps[0]);CHKERRQ(ierr);
6370f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps up           %d\n",jac->gridsweeps[1]);CHKERRQ(ierr);
6380f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on coarse    %d\n",jac->gridsweeps[2]);CHKERRQ(ierr);
63916d9e3a6SLisandro Dalcin 
6400f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr);
6410f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr);
6420f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr);
64316d9e3a6SLisandro Dalcin 
64416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax weight  (all)      %G\n",jac->relaxweight);CHKERRQ(ierr);
64516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr);
64616d9e3a6SLisandro Dalcin 
64716d9e3a6SLisandro Dalcin     if (jac->relaxorder) {
64816d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr);
64916d9e3a6SLisandro Dalcin     } else {
65016d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr);
65116d9e3a6SLisandro Dalcin     }
65216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr);
65316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr);
6540f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation type  %s\n",HYPREBoomerAMGInterpType[jac->interptype]);CHKERRQ(ierr);
6558f87f92bSBarry Smith     if (jac->nodal_coarsen) {
6568f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal())\n");CHKERRQ(ierr);
6578f87f92bSBarry Smith     }
6588f87f92bSBarry Smith     if (jac->nodal_relax) {
6598f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);CHKERRQ(ierr);
6608f87f92bSBarry Smith     }
66116d9e3a6SLisandro Dalcin   }
66216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
66316d9e3a6SLisandro Dalcin }
66416d9e3a6SLisandro Dalcin 
66516d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
66616d9e3a6SLisandro Dalcin #undef __FUNCT__
66716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails"
66816d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc)
66916d9e3a6SLisandro Dalcin {
67016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
67116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
6724ddd07fcSJed Brown   PetscInt       indx;
673ace3abfcSBarry Smith   PetscBool      flag;
67416d9e3a6SLisandro Dalcin   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};
67516d9e3a6SLisandro Dalcin 
67616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
67716d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr);
67816d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr);
67916d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr);
6802fa5cd67SKarl Rupp   if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels));
68116d9e3a6SLisandro Dalcin 
68216d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr);
6832fa5cd67SKarl Rupp   if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter));
68416d9e3a6SLisandro Dalcin 
68516d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr);
6862fa5cd67SKarl Rupp   if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal));
68716d9e3a6SLisandro Dalcin 
688acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_parasails_logging","Print info to screen","None",(PetscBool)jac->logging,(PetscBool*)&jac->logging,&flag);CHKERRQ(ierr);
6892fa5cd67SKarl Rupp   if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging));
69016d9e3a6SLisandro Dalcin 
691acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscBool)jac->ruse,(PetscBool*)&jac->ruse,&flag);CHKERRQ(ierr);
6922fa5cd67SKarl Rupp   if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse));
69316d9e3a6SLisandro Dalcin 
694a669f990SJed Brown   ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,ALEN(symtlist),symtlist[0],&indx,&flag);CHKERRQ(ierr);
69516d9e3a6SLisandro Dalcin   if (flag) {
69616d9e3a6SLisandro Dalcin     jac->symt = indx;
697fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
69816d9e3a6SLisandro Dalcin   }
69916d9e3a6SLisandro Dalcin 
70016d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
70116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
70216d9e3a6SLisandro Dalcin }
70316d9e3a6SLisandro Dalcin 
70416d9e3a6SLisandro Dalcin #undef __FUNCT__
70516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_ParaSails"
70616d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
70716d9e3a6SLisandro Dalcin {
70816d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
70916d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
710ace3abfcSBarry Smith   PetscBool      iascii;
71116d9e3a6SLisandro Dalcin   const char     *symt = 0;;
71216d9e3a6SLisandro Dalcin 
71316d9e3a6SLisandro Dalcin   PetscFunctionBegin;
714251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
71516d9e3a6SLisandro Dalcin   if (iascii) {
71616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");CHKERRQ(ierr);
71716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr);
71816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr);
71916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr);
72016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr);
721ace3abfcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscBools[jac->ruse]);CHKERRQ(ierr);
722ace3abfcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscBools[jac->logging]);CHKERRQ(ierr);
7232fa5cd67SKarl Rupp     if (!jac->symt) symt = "nonsymmetric matrix and preconditioner";
7242fa5cd67SKarl Rupp     else if (jac->symt == 1) symt = "SPD matrix and preconditioner";
7252fa5cd67SKarl Rupp     else if (jac->symt == 2) symt = "nonsymmetric matrix but SPD preconditioner";
726ce94432eSBarry Smith     else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
72716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr);
72816d9e3a6SLisandro Dalcin   }
72916d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
73016d9e3a6SLisandro Dalcin }
73116d9e3a6SLisandro Dalcin /* ---------------------------------------------------------------------------------*/
73216d9e3a6SLisandro Dalcin 
73316d9e3a6SLisandro Dalcin #undef __FUNCT__
73416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType_HYPRE"
735f7a08781SBarry Smith static PetscErrorCode  PCHYPREGetType_HYPRE(PC pc,const char *name[])
73616d9e3a6SLisandro Dalcin {
73716d9e3a6SLisandro Dalcin   PC_HYPRE *jac = (PC_HYPRE*)pc->data;
73816d9e3a6SLisandro Dalcin 
73916d9e3a6SLisandro Dalcin   PetscFunctionBegin;
74016d9e3a6SLisandro Dalcin   *name = jac->hypre_type;
74116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
74216d9e3a6SLisandro Dalcin }
74316d9e3a6SLisandro Dalcin 
74416d9e3a6SLisandro Dalcin #undef __FUNCT__
74516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType_HYPRE"
746f7a08781SBarry Smith static PetscErrorCode  PCHYPRESetType_HYPRE(PC pc,const char name[])
74716d9e3a6SLisandro Dalcin {
74816d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
74916d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
750ace3abfcSBarry Smith   PetscBool      flag;
75116d9e3a6SLisandro Dalcin 
75216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
75316d9e3a6SLisandro Dalcin   if (jac->hypre_type) {
75416d9e3a6SLisandro Dalcin     ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr);
755ce94432eSBarry Smith     if (!flag) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
75616d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
75716d9e3a6SLisandro Dalcin   } else {
75816d9e3a6SLisandro Dalcin     ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr);
75916d9e3a6SLisandro Dalcin   }
76016d9e3a6SLisandro Dalcin 
76116d9e3a6SLisandro Dalcin   jac->maxiter         = PETSC_DEFAULT;
76216d9e3a6SLisandro Dalcin   jac->tol             = PETSC_DEFAULT;
76316d9e3a6SLisandro Dalcin   jac->printstatistics = PetscLogPrintInfo;
76416d9e3a6SLisandro Dalcin 
76516d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr);
76616d9e3a6SLisandro Dalcin   if (flag) {
767fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParCSRPilutCreate,(jac->comm_hypre,&jac->hsolver));
76816d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
76916d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Pilut;
77016d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParCSRPilutDestroy;
77116d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParCSRPilutSetup;
77216d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParCSRPilutSolve;
77316d9e3a6SLisandro Dalcin     jac->factorrowsize      = PETSC_DEFAULT;
77416d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
77516d9e3a6SLisandro Dalcin   }
77616d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr);
77716d9e3a6SLisandro Dalcin   if (flag) {
778fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsCreate,(jac->comm_hypre,&jac->hsolver));
77916d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
78016d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_ParaSails;
78116d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParaSailsDestroy;
78216d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParaSailsSetup;
78316d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParaSailsSolve;
78416d9e3a6SLisandro Dalcin     /* initialize */
78516d9e3a6SLisandro Dalcin     jac->nlevels    = 1;
78616d9e3a6SLisandro Dalcin     jac->threshhold = .1;
78716d9e3a6SLisandro Dalcin     jac->filter     = .1;
78816d9e3a6SLisandro Dalcin     jac->loadbal    = 0;
7892fa5cd67SKarl Rupp     if (PetscLogPrintInfo) jac->logging = (int) PETSC_TRUE;
7902fa5cd67SKarl Rupp     else jac->logging = (int) PETSC_FALSE;
7912fa5cd67SKarl Rupp 
79216d9e3a6SLisandro Dalcin     jac->ruse = (int) PETSC_FALSE;
79316d9e3a6SLisandro Dalcin     jac->symt = 0;
794fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels));
795fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter));
796fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal));
797fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging));
798fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse));
799fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
80016d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
80116d9e3a6SLisandro Dalcin   }
80216d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr);
80316d9e3a6SLisandro Dalcin   if (flag) {
80416d9e3a6SLisandro Dalcin     ierr                    = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver);
80516d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
80616d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Euclid;
80716d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_EuclidDestroy;
80816d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_EuclidSetup;
80916d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_EuclidSolve;
81016d9e3a6SLisandro Dalcin     /* initialization */
81116d9e3a6SLisandro Dalcin     jac->bjilu              = PETSC_FALSE;
81216d9e3a6SLisandro Dalcin     jac->levels             = 1;
81316d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
81416d9e3a6SLisandro Dalcin   }
81516d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr);
81616d9e3a6SLisandro Dalcin   if (flag) {
81716d9e3a6SLisandro Dalcin     ierr                     = HYPRE_BoomerAMGCreate(&jac->hsolver);
81816d9e3a6SLisandro Dalcin     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
81916d9e3a6SLisandro Dalcin     pc->ops->view            = PCView_HYPRE_BoomerAMG;
82016d9e3a6SLisandro Dalcin     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
82116d9e3a6SLisandro Dalcin     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
82216d9e3a6SLisandro Dalcin     jac->destroy             = HYPRE_BoomerAMGDestroy;
82316d9e3a6SLisandro Dalcin     jac->setup               = HYPRE_BoomerAMGSetup;
82416d9e3a6SLisandro Dalcin     jac->solve               = HYPRE_BoomerAMGSolve;
82516d9e3a6SLisandro Dalcin     jac->applyrichardson     = PETSC_FALSE;
82616d9e3a6SLisandro Dalcin     /* these defaults match the hypre defaults */
82716d9e3a6SLisandro Dalcin     jac->cycletype        = 1;
82816d9e3a6SLisandro Dalcin     jac->maxlevels        = 25;
82916d9e3a6SLisandro Dalcin     jac->maxiter          = 1;
8308f87f92bSBarry Smith     jac->tol              = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */
83116d9e3a6SLisandro Dalcin     jac->truncfactor      = 0.0;
83216d9e3a6SLisandro Dalcin     jac->strongthreshold  = .25;
83316d9e3a6SLisandro Dalcin     jac->maxrowsum        = .9;
83416d9e3a6SLisandro Dalcin     jac->coarsentype      = 6;
83516d9e3a6SLisandro Dalcin     jac->measuretype      = 0;
8360f1074feSSatish Balay     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = 1;
8378f87f92bSBarry Smith     jac->relaxtype[0]     = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */
8380f1074feSSatish Balay     jac->relaxtype[2]     = 9; /*G.E. */
83916d9e3a6SLisandro Dalcin     jac->relaxweight      = 1.0;
84016d9e3a6SLisandro Dalcin     jac->outerrelaxweight = 1.0;
84116d9e3a6SLisandro Dalcin     jac->relaxorder       = 1;
8420f1074feSSatish Balay     jac->interptype       = 0;
8430f1074feSSatish Balay     jac->agg_nl           = 0;
8440f1074feSSatish Balay     jac->pmax             = 0;
8450f1074feSSatish Balay     jac->truncfactor      = 0.0;
8460f1074feSSatish Balay     jac->agg_num_paths    = 1;
8478f87f92bSBarry Smith 
8488f87f92bSBarry Smith     jac->nodal_coarsen      = 0;
8498f87f92bSBarry Smith     jac->nodal_relax        = PETSC_FALSE;
8508f87f92bSBarry Smith     jac->nodal_relax_levels = 1;
851fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
852fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
853fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
854fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
855fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
856fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
857fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
858fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
859fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
860fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
861fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
862fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
863fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
864fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
865fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, jac->relaxtype[0]));  /*defaults coarse to 9*/
866fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver, jac->gridsweeps[0])); /*defaults coarse to 1 */
86716d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
86816d9e3a6SLisandro Dalcin   }
869503cfb0cSBarry Smith   ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr);
8702fa5cd67SKarl Rupp 
8710298fd71SBarry Smith   jac->hypre_type = NULL;
872ce94432eSBarry Smith   SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
87316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
87416d9e3a6SLisandro Dalcin }
87516d9e3a6SLisandro Dalcin 
87616d9e3a6SLisandro Dalcin /*
87716d9e3a6SLisandro Dalcin     It only gets here if the HYPRE type has not been set before the call to
87816d9e3a6SLisandro Dalcin    ...SetFromOptions() which actually is most of the time
87916d9e3a6SLisandro Dalcin */
88016d9e3a6SLisandro Dalcin #undef __FUNCT__
88116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE"
88216d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE(PC pc)
88316d9e3a6SLisandro Dalcin {
88416d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
8854ddd07fcSJed Brown   PetscInt       indx;
88616d9e3a6SLisandro Dalcin   const char     *type[] = {"pilut","parasails","boomeramg","euclid"};
887ace3abfcSBarry Smith   PetscBool      flg;
88816d9e3a6SLisandro Dalcin 
88916d9e3a6SLisandro Dalcin   PetscFunctionBegin;
89016d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr);
89102a17cd4SBarry Smith   ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"boomeramg",&indx,&flg);CHKERRQ(ierr);
89216d9e3a6SLisandro Dalcin   if (flg) {
89316d9e3a6SLisandro Dalcin     ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr);
89402a17cd4SBarry Smith   } else {
89502a17cd4SBarry Smith     ierr = PCHYPRESetType_HYPRE(pc,"boomeramg");CHKERRQ(ierr);
89616d9e3a6SLisandro Dalcin   }
89716d9e3a6SLisandro Dalcin   if (pc->ops->setfromoptions) {
89816d9e3a6SLisandro Dalcin     ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr);
89916d9e3a6SLisandro Dalcin   }
90016d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
90116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
90216d9e3a6SLisandro Dalcin }
90316d9e3a6SLisandro Dalcin 
90416d9e3a6SLisandro Dalcin #undef __FUNCT__
90516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType"
90616d9e3a6SLisandro Dalcin /*@C
90716d9e3a6SLisandro Dalcin      PCHYPRESetType - Sets which hypre preconditioner you wish to use
90816d9e3a6SLisandro Dalcin 
90916d9e3a6SLisandro Dalcin    Input Parameters:
91016d9e3a6SLisandro Dalcin +     pc - the preconditioner context
91116d9e3a6SLisandro Dalcin -     name - either  pilut, parasails, boomeramg, euclid
91216d9e3a6SLisandro Dalcin 
91316d9e3a6SLisandro Dalcin    Options Database Keys:
91416d9e3a6SLisandro Dalcin    -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
91516d9e3a6SLisandro Dalcin 
91616d9e3a6SLisandro Dalcin    Level: intermediate
91716d9e3a6SLisandro Dalcin 
91816d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
91916d9e3a6SLisandro Dalcin            PCHYPRE
92016d9e3a6SLisandro Dalcin 
92116d9e3a6SLisandro Dalcin @*/
9227087cfbeSBarry Smith PetscErrorCode  PCHYPRESetType(PC pc,const char name[])
92316d9e3a6SLisandro Dalcin {
9244ac538c5SBarry Smith   PetscErrorCode ierr;
92516d9e3a6SLisandro Dalcin 
92616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9270700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
92816d9e3a6SLisandro Dalcin   PetscValidCharPointer(name,2);
9294ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCHYPRESetType_C",(PC,const char[]),(pc,name));CHKERRQ(ierr);
93016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
93116d9e3a6SLisandro Dalcin }
93216d9e3a6SLisandro Dalcin 
93316d9e3a6SLisandro Dalcin #undef __FUNCT__
93416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType"
93516d9e3a6SLisandro Dalcin /*@C
93616d9e3a6SLisandro Dalcin      PCHYPREGetType - Gets which hypre preconditioner you are using
93716d9e3a6SLisandro Dalcin 
93816d9e3a6SLisandro Dalcin    Input Parameter:
93916d9e3a6SLisandro Dalcin .     pc - the preconditioner context
94016d9e3a6SLisandro Dalcin 
94116d9e3a6SLisandro Dalcin    Output Parameter:
94216d9e3a6SLisandro Dalcin .     name - either  pilut, parasails, boomeramg, euclid
94316d9e3a6SLisandro Dalcin 
94416d9e3a6SLisandro Dalcin    Level: intermediate
94516d9e3a6SLisandro Dalcin 
94616d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC,
94716d9e3a6SLisandro Dalcin            PCHYPRE
94816d9e3a6SLisandro Dalcin 
94916d9e3a6SLisandro Dalcin @*/
9507087cfbeSBarry Smith PetscErrorCode  PCHYPREGetType(PC pc,const char *name[])
95116d9e3a6SLisandro Dalcin {
9524ac538c5SBarry Smith   PetscErrorCode ierr;
95316d9e3a6SLisandro Dalcin 
95416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9550700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
95616d9e3a6SLisandro Dalcin   PetscValidPointer(name,2);
9574ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCHYPREGetType_C",(PC,const char*[]),(pc,name));CHKERRQ(ierr);
95816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
95916d9e3a6SLisandro Dalcin }
96016d9e3a6SLisandro Dalcin 
96116d9e3a6SLisandro Dalcin /*MC
96216d9e3a6SLisandro Dalcin      PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre
96316d9e3a6SLisandro Dalcin 
96416d9e3a6SLisandro Dalcin    Options Database Keys:
96516d9e3a6SLisandro Dalcin +   -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
96616d9e3a6SLisandro Dalcin -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
96716d9e3a6SLisandro Dalcin           preconditioner
96816d9e3a6SLisandro Dalcin 
96916d9e3a6SLisandro Dalcin    Level: intermediate
97016d9e3a6SLisandro Dalcin 
97116d9e3a6SLisandro Dalcin    Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()),
97216d9e3a6SLisandro Dalcin           the many hypre options can ONLY be set via the options database (e.g. the command line
97316d9e3a6SLisandro Dalcin           or with PetscOptionsSetValue(), there are no functions to set them)
97416d9e3a6SLisandro Dalcin 
97516d9e3a6SLisandro Dalcin           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations
9760f1074feSSatish Balay           (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if
9770f1074feSSatish Balay           -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner
9780f1074feSSatish Balay           (-pc_hypre_boomeramg_rtol should be set to 0.0 - the default - to strictly use a fixed number of
9798f87f92bSBarry Smith           iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations
9800f1074feSSatish Balay           and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10
9810f1074feSSatish Balay           then AT MOST twenty V-cycles of boomeramg will be called.
98216d9e3a6SLisandro Dalcin 
9830f1074feSSatish Balay            Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation
9840f1074feSSatish Balay            (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry.
9850f1074feSSatish Balay            Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi.
98616d9e3a6SLisandro Dalcin           If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly
98716d9e3a6SLisandro Dalcin           and use -ksp_max_it to control the number of V-cycles.
98816d9e3a6SLisandro Dalcin           (see the PETSc FAQ.html at the PETSc website under the Documentation tab).
98916d9e3a6SLisandro Dalcin 
99016d9e3a6SLisandro Dalcin           2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option
99116d9e3a6SLisandro Dalcin           -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L.
99216d9e3a6SLisandro Dalcin 
9939e5bc791SBarry Smith           See PCPFMG for access to the hypre Struct PFMG solver
9949e5bc791SBarry Smith 
99516d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
9969e5bc791SBarry Smith            PCHYPRESetType(), PCPFMG
99716d9e3a6SLisandro Dalcin 
99816d9e3a6SLisandro Dalcin M*/
99916d9e3a6SLisandro Dalcin 
100016d9e3a6SLisandro Dalcin #undef __FUNCT__
100116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCCreate_HYPRE"
10028cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_HYPRE(PC pc)
100316d9e3a6SLisandro Dalcin {
100416d9e3a6SLisandro Dalcin   PC_HYPRE       *jac;
100516d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
100616d9e3a6SLisandro Dalcin 
100716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
1008*b00a9115SJed Brown   ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr);
10092fa5cd67SKarl Rupp 
101016d9e3a6SLisandro Dalcin   pc->data                = jac;
101116d9e3a6SLisandro Dalcin   pc->ops->destroy        = PCDestroy_HYPRE;
101216d9e3a6SLisandro Dalcin   pc->ops->setfromoptions = PCSetFromOptions_HYPRE;
101316d9e3a6SLisandro Dalcin   pc->ops->setup          = PCSetUp_HYPRE;
101416d9e3a6SLisandro Dalcin   pc->ops->apply          = PCApply_HYPRE;
101516d9e3a6SLisandro Dalcin   jac->comm_hypre         = MPI_COMM_NULL;
10160298fd71SBarry Smith   jac->hypre_type         = NULL;
101716d9e3a6SLisandro Dalcin   /* duplicate communicator for hypre */
1018ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(jac->comm_hypre));CHKERRQ(ierr);
1019bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetType_C",PCHYPRESetType_HYPRE);CHKERRQ(ierr);
1020bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHYPREGetType_C",PCHYPREGetType_HYPRE);CHKERRQ(ierr);
102116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
102216d9e3a6SLisandro Dalcin }
1023ebc551c0SBarry Smith 
1024f91d8e95SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------*/
1025f91d8e95SBarry Smith 
1026b862ddfaSBarry Smith /* this include is needed ONLY to allow access to the private data inside the Mat object specific to hypre */
1027b45d2f2cSJed Brown #include <petsc-private/matimpl.h>
1028ebc551c0SBarry Smith 
1029ebc551c0SBarry Smith typedef struct {
103068326731SBarry Smith   MPI_Comm           hcomm;        /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */
1031f91d8e95SBarry Smith   HYPRE_StructSolver hsolver;
10329e5bc791SBarry Smith 
10339e5bc791SBarry Smith   /* keep copy of PFMG options used so may view them */
10344ddd07fcSJed Brown   PetscInt its;
10359e5bc791SBarry Smith   double   tol;
10364ddd07fcSJed Brown   PetscInt relax_type;
10374ddd07fcSJed Brown   PetscInt rap_type;
10384ddd07fcSJed Brown   PetscInt num_pre_relax,num_post_relax;
10394ddd07fcSJed Brown   PetscInt max_levels;
1040ebc551c0SBarry Smith } PC_PFMG;
1041ebc551c0SBarry Smith 
1042ebc551c0SBarry Smith #undef __FUNCT__
1043ebc551c0SBarry Smith #define __FUNCT__ "PCDestroy_PFMG"
1044ebc551c0SBarry Smith PetscErrorCode PCDestroy_PFMG(PC pc)
1045ebc551c0SBarry Smith {
1046ebc551c0SBarry Smith   PetscErrorCode ierr;
1047f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1048ebc551c0SBarry Smith 
1049ebc551c0SBarry Smith   PetscFunctionBegin;
10502fa5cd67SKarl Rupp   if (ex->hsolver) PetscStackCallStandard(HYPRE_StructPFMGDestroy,(ex->hsolver));
1051f91d8e95SBarry Smith   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1052c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1053ebc551c0SBarry Smith   PetscFunctionReturn(0);
1054ebc551c0SBarry Smith }
1055ebc551c0SBarry Smith 
10569e5bc791SBarry Smith static const char *PFMGRelaxType[] = {"Jacobi","Weighted-Jacobi","symmetric-Red/Black-Gauss-Seidel","Red/Black-Gauss-Seidel"};
10579e5bc791SBarry Smith static const char *PFMGRAPType[] = {"Galerkin","non-Galerkin"};
10589e5bc791SBarry Smith 
1059ebc551c0SBarry Smith #undef __FUNCT__
1060ebc551c0SBarry Smith #define __FUNCT__ "PCView_PFMG"
1061ebc551c0SBarry Smith PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer)
1062ebc551c0SBarry Smith {
1063ebc551c0SBarry Smith   PetscErrorCode ierr;
1064ace3abfcSBarry Smith   PetscBool      iascii;
1065f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1066ebc551c0SBarry Smith 
1067ebc551c0SBarry Smith   PetscFunctionBegin;
1068251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10699e5bc791SBarry Smith   if (iascii) {
10709e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG preconditioning\n");CHKERRQ(ierr);
10719e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
10729e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
10739e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
10749e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: RAP type %s\n",PFMGRAPType[ex->rap_type]);CHKERRQ(ierr);
10759e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
10763b46a515SGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max levels %d\n",ex->max_levels);CHKERRQ(ierr);
10779e5bc791SBarry Smith   }
1078ebc551c0SBarry Smith   PetscFunctionReturn(0);
1079ebc551c0SBarry Smith }
1080ebc551c0SBarry Smith 
10819e5bc791SBarry Smith 
1082ebc551c0SBarry Smith #undef __FUNCT__
1083ebc551c0SBarry Smith #define __FUNCT__ "PCSetFromOptions_PFMG"
1084ebc551c0SBarry Smith PetscErrorCode PCSetFromOptions_PFMG(PC pc)
1085ebc551c0SBarry Smith {
1086ebc551c0SBarry Smith   PetscErrorCode ierr;
1087f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1088ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
1089ebc551c0SBarry Smith 
1090ebc551c0SBarry Smith   PetscFunctionBegin;
1091ebc551c0SBarry Smith   ierr = PetscOptionsHead("PFMG options");CHKERRQ(ierr);
10920298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",flg,&flg,NULL);CHKERRQ(ierr);
109368326731SBarry Smith   if (flg) {
1094a0324ebeSBarry Smith     int level=3;
1095fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_StructPFMGSetPrintLevel,(ex->hsolver,level));
109668326731SBarry Smith   }
10970298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_pfmg_its","Number of iterations of PFMG to use as preconditioner","HYPRE_StructPFMGSetMaxIter",ex->its,&ex->its,NULL);CHKERRQ(ierr);
1098fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(ex->hsolver,ex->its));
10990298fd71SBarry 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,NULL);CHKERRQ(ierr);
1100fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetNumPreRelax,(ex->hsolver,ex->num_pre_relax));
11010298fd71SBarry 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,NULL);CHKERRQ(ierr);
1102fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetNumPostRelax,(ex->hsolver,ex->num_post_relax));
11039e5bc791SBarry Smith 
11040298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_pfmg_max_levels","Max Levels for MG hierarchy","HYPRE_StructPFMGSetMaxLevels",ex->max_levels,&ex->max_levels,NULL);CHKERRQ(ierr);
1105fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetMaxLevels,(ex->hsolver,ex->max_levels));
11063b46a515SGlenn Hammond 
11070298fd71SBarry Smith   ierr = PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,NULL);CHKERRQ(ierr);
1108fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(ex->hsolver,ex->tol));
11090298fd71SBarry Smith   ierr = PetscOptionsEList("-pc_pfmg_relax_type","Relax type for the up and down cycles","HYPRE_StructPFMGSetRelaxType",PFMGRelaxType,ALEN(PFMGRelaxType),PFMGRelaxType[ex->relax_type],&ex->relax_type,NULL);CHKERRQ(ierr);
1110fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetRelaxType,(ex->hsolver, ex->relax_type));
11110298fd71SBarry Smith   ierr = PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,ALEN(PFMGRAPType),PFMGRAPType[ex->rap_type],&ex->rap_type,NULL);CHKERRQ(ierr);
1112fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetRAPType,(ex->hsolver, ex->rap_type));
1113ebc551c0SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1114ebc551c0SBarry Smith   PetscFunctionReturn(0);
1115ebc551c0SBarry Smith }
1116ebc551c0SBarry Smith 
1117f91d8e95SBarry Smith #undef __FUNCT__
1118f91d8e95SBarry Smith #define __FUNCT__ "PCApply_PFMG"
1119f91d8e95SBarry Smith PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y)
1120f91d8e95SBarry Smith {
1121f91d8e95SBarry Smith   PetscErrorCode  ierr;
1122f91d8e95SBarry Smith   PC_PFMG         *ex = (PC_PFMG*) pc->data;
1123f91d8e95SBarry Smith   PetscScalar     *xx,*yy;
11244ddd07fcSJed Brown   PetscInt        ilower[3],iupper[3];
112568326731SBarry Smith   Mat_HYPREStruct *mx = (Mat_HYPREStruct*)(pc->pmat->data);
1126f91d8e95SBarry Smith 
1127f91d8e95SBarry Smith   PetscFunctionBegin;
1128dff31646SBarry Smith   ierr       = PetscCitationsRegister(hypreCitation,&cite);CHKERRQ(ierr);
1129aa219208SBarry Smith   ierr       = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1130f91d8e95SBarry Smith   iupper[0] += ilower[0] - 1;
1131f91d8e95SBarry Smith   iupper[1] += ilower[1] - 1;
1132f91d8e95SBarry Smith   iupper[2] += ilower[2] - 1;
1133f91d8e95SBarry Smith 
1134f91d8e95SBarry Smith   /* copy x values over to hypre */
1135fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorSetConstantValues,(mx->hb,0.0));
1136f91d8e95SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
11378b1f7689SBarry Smith   PetscStackCallStandard(HYPRE_StructVectorSetBoxValues,(mx->hb,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,xx));
1138f91d8e95SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1139fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorAssemble,(mx->hb));
1140fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSolve,(ex->hsolver,mx->hmat,mx->hb,mx->hx));
1141f91d8e95SBarry Smith 
1142f91d8e95SBarry Smith   /* copy solution values back to PETSc */
1143f91d8e95SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
11448b1f7689SBarry Smith   PetscStackCallStandard(HYPRE_StructVectorGetBoxValues,(mx->hx,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,yy));
1145f91d8e95SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1146f91d8e95SBarry Smith   PetscFunctionReturn(0);
1147f91d8e95SBarry Smith }
1148f91d8e95SBarry Smith 
11499e5bc791SBarry Smith #undef __FUNCT__
11509e5bc791SBarry Smith #define __FUNCT__ "PCApplyRichardson_PFMG"
1151ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_PFMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
11529e5bc791SBarry Smith {
11539e5bc791SBarry Smith   PC_PFMG        *jac = (PC_PFMG*)pc->data;
11549e5bc791SBarry Smith   PetscErrorCode ierr;
11554ddd07fcSJed Brown   PetscInt       oits;
11569e5bc791SBarry Smith 
11579e5bc791SBarry Smith   PetscFunctionBegin;
1158dff31646SBarry Smith   ierr = PetscCitationsRegister(hypreCitation,&cite);CHKERRQ(ierr);
1159fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(jac->hsolver,its*jac->its));
1160fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(jac->hsolver,rtol));
11619e5bc791SBarry Smith 
11629e5bc791SBarry Smith   ierr = PCApply_PFMG(pc,b,y);CHKERRQ(ierr);
11638b1f7689SBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGGetNumIterations,(jac->hsolver,(HYPRE_Int *)&oits));
11649e5bc791SBarry Smith   *outits = oits;
11659e5bc791SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
11669e5bc791SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
1167fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetTol,(jac->hsolver,jac->tol));
1168fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(jac->hsolver,jac->its));
11699e5bc791SBarry Smith   PetscFunctionReturn(0);
11709e5bc791SBarry Smith }
11719e5bc791SBarry Smith 
11729e5bc791SBarry Smith 
11733a32d3dbSGlenn Hammond #undef __FUNCT__
11743a32d3dbSGlenn Hammond #define __FUNCT__ "PCSetUp_PFMG"
11753a32d3dbSGlenn Hammond PetscErrorCode PCSetUp_PFMG(PC pc)
11763a32d3dbSGlenn Hammond {
11773a32d3dbSGlenn Hammond   PetscErrorCode  ierr;
11783a32d3dbSGlenn Hammond   PC_PFMG         *ex = (PC_PFMG*) pc->data;
11793a32d3dbSGlenn Hammond   Mat_HYPREStruct *mx = (Mat_HYPREStruct*)(pc->pmat->data);
1180ace3abfcSBarry Smith   PetscBool       flg;
11813a32d3dbSGlenn Hammond 
11823a32d3dbSGlenn Hammond   PetscFunctionBegin;
1183251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESTRUCT,&flg);CHKERRQ(ierr);
1184ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESTRUCT with this preconditioner");
11853a32d3dbSGlenn Hammond 
11863a32d3dbSGlenn Hammond   /* create the hypre solver object and set its information */
11872fa5cd67SKarl Rupp   if (ex->hsolver) PetscStackCallStandard(HYPRE_StructPFMGDestroy,(ex->hsolver));
1188fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
11893a32d3dbSGlenn Hammond   ierr = PCSetFromOptions_PFMG(pc);CHKERRQ(ierr);
1190fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetup,(ex->hsolver,mx->hmat,mx->hb,mx->hx));
1191fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGSetZeroGuess,(ex->hsolver));
11923a32d3dbSGlenn Hammond   PetscFunctionReturn(0);
11933a32d3dbSGlenn Hammond }
11943a32d3dbSGlenn Hammond 
1195ebc551c0SBarry Smith 
1196ebc551c0SBarry Smith /*MC
1197ebc551c0SBarry Smith      PCPFMG - the hypre PFMG multigrid solver
1198ebc551c0SBarry Smith 
1199ebc551c0SBarry Smith    Level: advanced
1200ebc551c0SBarry Smith 
12019e5bc791SBarry Smith    Options Database:
12029e5bc791SBarry Smith + -pc_pfmg_its <its> number of iterations of PFMG to use as preconditioner
12039e5bc791SBarry Smith . -pc_pfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
12049e5bc791SBarry Smith . -pc_pfmg_num_post_relax <steps> number of smoothing steps after coarse grid
12059e5bc791SBarry Smith . -pc_pfmg_tol <tol> tolerance of PFMG
12069e5bc791SBarry 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
12079e5bc791SBarry Smith - -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin
1208f91d8e95SBarry Smith 
12099e5bc791SBarry Smith    Notes:  This is for CELL-centered descretizations
12109e5bc791SBarry Smith 
12118e395302SJed Brown            This must be used with the MATHYPRESTRUCT matrix type.
1212aa219208SBarry Smith            This is less general than in hypre, it supports only one block per process defined by a PETSc DMDA.
12139e5bc791SBarry Smith 
12149e5bc791SBarry Smith .seealso:  PCMG, MATHYPRESTRUCT
1215ebc551c0SBarry Smith M*/
1216ebc551c0SBarry Smith 
1217ebc551c0SBarry Smith #undef __FUNCT__
1218ebc551c0SBarry Smith #define __FUNCT__ "PCCreate_PFMG"
12198cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_PFMG(PC pc)
1220ebc551c0SBarry Smith {
1221ebc551c0SBarry Smith   PetscErrorCode ierr;
1222ebc551c0SBarry Smith   PC_PFMG        *ex;
1223ebc551c0SBarry Smith 
1224ebc551c0SBarry Smith   PetscFunctionBegin;
1225*b00a9115SJed Brown   ierr     = PetscNew(&ex);CHKERRQ(ierr); \
122668326731SBarry Smith   pc->data = ex;
1227ebc551c0SBarry Smith 
12289e5bc791SBarry Smith   ex->its            = 1;
12299e5bc791SBarry Smith   ex->tol            = 1.e-8;
12309e5bc791SBarry Smith   ex->relax_type     = 1;
12319e5bc791SBarry Smith   ex->rap_type       = 0;
12329e5bc791SBarry Smith   ex->num_pre_relax  = 1;
12339e5bc791SBarry Smith   ex->num_post_relax = 1;
12343b46a515SGlenn Hammond   ex->max_levels     = 0;
12359e5bc791SBarry Smith 
1236ebc551c0SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_PFMG;
1237ebc551c0SBarry Smith   pc->ops->view            = PCView_PFMG;
1238ebc551c0SBarry Smith   pc->ops->destroy         = PCDestroy_PFMG;
1239f91d8e95SBarry Smith   pc->ops->apply           = PCApply_PFMG;
12409e5bc791SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_PFMG;
124168326731SBarry Smith   pc->ops->setup           = PCSetUp_PFMG;
12422fa5cd67SKarl Rupp 
1243ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ex->hcomm));CHKERRQ(ierr);
1244fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
1245ebc551c0SBarry Smith   PetscFunctionReturn(0);
1246ebc551c0SBarry Smith }
1247d851a50bSGlenn Hammond 
1248325fc9f4SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------------------------*/
1249325fc9f4SBarry Smith 
1250d851a50bSGlenn Hammond /* we know we are working with a HYPRE_SStructMatrix */
1251d851a50bSGlenn Hammond typedef struct {
1252d851a50bSGlenn Hammond   MPI_Comm            hcomm;       /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */
1253d851a50bSGlenn Hammond   HYPRE_SStructSolver ss_solver;
1254d851a50bSGlenn Hammond 
1255d851a50bSGlenn Hammond   /* keep copy of SYSPFMG options used so may view them */
12564ddd07fcSJed Brown   PetscInt its;
1257d851a50bSGlenn Hammond   double   tol;
12584ddd07fcSJed Brown   PetscInt relax_type;
12594ddd07fcSJed Brown   PetscInt num_pre_relax,num_post_relax;
1260d851a50bSGlenn Hammond } PC_SysPFMG;
1261d851a50bSGlenn Hammond 
1262d851a50bSGlenn Hammond #undef __FUNCT__
1263d851a50bSGlenn Hammond #define __FUNCT__ "PCDestroy_SysPFMG"
1264d851a50bSGlenn Hammond PetscErrorCode PCDestroy_SysPFMG(PC pc)
1265d851a50bSGlenn Hammond {
1266d851a50bSGlenn Hammond   PetscErrorCode ierr;
1267d851a50bSGlenn Hammond   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;
1268d851a50bSGlenn Hammond 
1269d851a50bSGlenn Hammond   PetscFunctionBegin;
12702fa5cd67SKarl Rupp   if (ex->ss_solver) PetscStackCallStandard(HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));
1271d851a50bSGlenn Hammond   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1272c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1273d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1274d851a50bSGlenn Hammond }
1275d851a50bSGlenn Hammond 
1276d851a50bSGlenn Hammond static const char *SysPFMGRelaxType[] = {"Weighted-Jacobi","Red/Black-Gauss-Seidel"};
1277d851a50bSGlenn Hammond 
1278d851a50bSGlenn Hammond #undef __FUNCT__
1279d851a50bSGlenn Hammond #define __FUNCT__ "PCView_SysPFMG"
1280d851a50bSGlenn Hammond PetscErrorCode PCView_SysPFMG(PC pc,PetscViewer viewer)
1281d851a50bSGlenn Hammond {
1282d851a50bSGlenn Hammond   PetscErrorCode ierr;
1283ace3abfcSBarry Smith   PetscBool      iascii;
1284d851a50bSGlenn Hammond   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;
1285d851a50bSGlenn Hammond 
1286d851a50bSGlenn Hammond   PetscFunctionBegin;
1287251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1288d851a50bSGlenn Hammond   if (iascii) {
1289d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG preconditioning\n");CHKERRQ(ierr);
1290d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
1291d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
1292d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
1293d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
1294d851a50bSGlenn Hammond   }
1295d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1296d851a50bSGlenn Hammond }
1297d851a50bSGlenn Hammond 
1298d851a50bSGlenn Hammond 
1299d851a50bSGlenn Hammond #undef __FUNCT__
1300d851a50bSGlenn Hammond #define __FUNCT__ "PCSetFromOptions_SysPFMG"
1301d851a50bSGlenn Hammond PetscErrorCode PCSetFromOptions_SysPFMG(PC pc)
1302d851a50bSGlenn Hammond {
1303d851a50bSGlenn Hammond   PetscErrorCode ierr;
1304d851a50bSGlenn Hammond   PC_SysPFMG     *ex = (PC_SysPFMG*) pc->data;
1305ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
1306d851a50bSGlenn Hammond 
1307d851a50bSGlenn Hammond   PetscFunctionBegin;
1308d851a50bSGlenn Hammond   ierr = PetscOptionsHead("SysPFMG options");CHKERRQ(ierr);
13090298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_syspfmg_print_statistics","Print statistics","HYPRE_SStructSysPFMGSetPrintLevel",flg,&flg,NULL);CHKERRQ(ierr);
1310d851a50bSGlenn Hammond   if (flg) {
1311d851a50bSGlenn Hammond     int level=3;
1312fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructSysPFMGSetPrintLevel,(ex->ss_solver,level));
1313d851a50bSGlenn Hammond   }
13140298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_syspfmg_its","Number of iterations of SysPFMG to use as preconditioner","HYPRE_SStructSysPFMGSetMaxIter",ex->its,&ex->its,NULL);CHKERRQ(ierr);
1315fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(ex->ss_solver,ex->its));
13160298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_syspfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_SStructSysPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,NULL);CHKERRQ(ierr);
1317fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetNumPreRelax,(ex->ss_solver,ex->num_pre_relax));
13180298fd71SBarry Smith   ierr = PetscOptionsInt("-pc_syspfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_SStructSysPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,NULL);CHKERRQ(ierr);
1319fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetNumPostRelax,(ex->ss_solver,ex->num_post_relax));
1320d851a50bSGlenn Hammond 
13210298fd71SBarry Smith   ierr = PetscOptionsReal("-pc_syspfmg_tol","Tolerance of SysPFMG","HYPRE_SStructSysPFMGSetTol",ex->tol,&ex->tol,NULL);CHKERRQ(ierr);
1322fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(ex->ss_solver,ex->tol));
13230298fd71SBarry Smith   ierr = PetscOptionsEList("-pc_syspfmg_relax_type","Relax type for the up and down cycles","HYPRE_SStructSysPFMGSetRelaxType",SysPFMGRelaxType,4,SysPFMGRelaxType[ex->relax_type],&ex->relax_type,NULL);CHKERRQ(ierr);
1324fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetRelaxType,(ex->ss_solver, ex->relax_type));
1325d851a50bSGlenn Hammond   ierr = PetscOptionsTail();CHKERRQ(ierr);
1326d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1327d851a50bSGlenn Hammond }
1328d851a50bSGlenn Hammond 
1329d851a50bSGlenn Hammond #undef __FUNCT__
1330d851a50bSGlenn Hammond #define __FUNCT__ "PCApply_SysPFMG"
1331d851a50bSGlenn Hammond PetscErrorCode PCApply_SysPFMG(PC pc,Vec x,Vec y)
1332d851a50bSGlenn Hammond {
1333d851a50bSGlenn Hammond   PetscErrorCode   ierr;
1334d851a50bSGlenn Hammond   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
1335d851a50bSGlenn Hammond   PetscScalar      *xx,*yy;
13364ddd07fcSJed Brown   PetscInt         ilower[3],iupper[3];
1337d851a50bSGlenn Hammond   Mat_HYPRESStruct *mx     = (Mat_HYPRESStruct*)(pc->pmat->data);
13384ddd07fcSJed Brown   PetscInt         ordering= mx->dofs_order;
13394ddd07fcSJed Brown   PetscInt         nvars   = mx->nvars;
13404ddd07fcSJed Brown   PetscInt         part    = 0;
13414ddd07fcSJed Brown   PetscInt         size;
13424ddd07fcSJed Brown   PetscInt         i;
1343d851a50bSGlenn Hammond 
1344d851a50bSGlenn Hammond   PetscFunctionBegin;
1345dff31646SBarry Smith   ierr       = PetscCitationsRegister(hypreCitation,&cite);CHKERRQ(ierr);
1346aa219208SBarry Smith   ierr       = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1347d851a50bSGlenn Hammond   iupper[0] += ilower[0] - 1;
1348d851a50bSGlenn Hammond   iupper[1] += ilower[1] - 1;
1349d851a50bSGlenn Hammond   iupper[2] += ilower[2] - 1;
1350d851a50bSGlenn Hammond 
1351d851a50bSGlenn Hammond   size = 1;
13522fa5cd67SKarl Rupp   for (i= 0; i< 3; i++) size *= (iupper[i]-ilower[i]+1);
13532fa5cd67SKarl Rupp 
1354d851a50bSGlenn Hammond   /* copy x values over to hypre for variable ordering */
1355d851a50bSGlenn Hammond   if (ordering) {
1356fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
1357d851a50bSGlenn Hammond     ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
13588b1f7689SBarry Smith     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,xx+(size*i)));
1359d851a50bSGlenn Hammond     ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1360fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
1361fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
1362fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1363d851a50bSGlenn Hammond 
1364d851a50bSGlenn Hammond     /* copy solution values back to PETSc */
1365d851a50bSGlenn Hammond     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
13668b1f7689SBarry Smith     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,yy+(size*i)));
1367d851a50bSGlenn Hammond     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1368a65764d7SBarry Smith   } else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
1369d851a50bSGlenn Hammond     PetscScalar *z;
13704ddd07fcSJed Brown     PetscInt    j, k;
1371d851a50bSGlenn Hammond 
1372785e854fSJed Brown     ierr = PetscMalloc1(nvars*size,&z);CHKERRQ(ierr);
1373fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
1374d851a50bSGlenn Hammond     ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1375d851a50bSGlenn Hammond 
1376d851a50bSGlenn Hammond     /* transform nodal to hypre's variable ordering for sys_pfmg */
1377d851a50bSGlenn Hammond     for (i= 0; i< size; i++) {
1378d851a50bSGlenn Hammond       k= i*nvars;
13792fa5cd67SKarl Rupp       for (j= 0; j< nvars; j++) z[j*size+i]= xx[k+j];
1380d851a50bSGlenn Hammond     }
13818b1f7689SBarry Smith     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,z+(size*i)));
1382d851a50bSGlenn Hammond     ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1383fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
1384fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1385d851a50bSGlenn Hammond 
1386d851a50bSGlenn Hammond     /* copy solution values back to PETSc */
1387d851a50bSGlenn Hammond     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
13888b1f7689SBarry Smith     for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,z+(size*i)));
1389d851a50bSGlenn Hammond     /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
1390d851a50bSGlenn Hammond     for (i= 0; i< size; i++) {
1391d851a50bSGlenn Hammond       k= i*nvars;
13922fa5cd67SKarl Rupp       for (j= 0; j< nvars; j++) yy[k+j]= z[j*size+i];
1393d851a50bSGlenn Hammond     }
1394d851a50bSGlenn Hammond     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1395d851a50bSGlenn Hammond     ierr = PetscFree(z);CHKERRQ(ierr);
1396d851a50bSGlenn Hammond   }
1397d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1398d851a50bSGlenn Hammond }
1399d851a50bSGlenn Hammond 
1400d851a50bSGlenn Hammond #undef __FUNCT__
1401d851a50bSGlenn Hammond #define __FUNCT__ "PCApplyRichardson_SysPFMG"
1402ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_SysPFMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
1403d851a50bSGlenn Hammond {
1404d851a50bSGlenn Hammond   PC_SysPFMG     *jac = (PC_SysPFMG*)pc->data;
1405d851a50bSGlenn Hammond   PetscErrorCode ierr;
14064ddd07fcSJed Brown   PetscInt       oits;
1407d851a50bSGlenn Hammond 
1408d851a50bSGlenn Hammond   PetscFunctionBegin;
1409dff31646SBarry Smith   ierr = PetscCitationsRegister(hypreCitation,&cite);CHKERRQ(ierr);
1410fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,its*jac->its));
1411fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,rtol));
1412d851a50bSGlenn Hammond   ierr = PCApply_SysPFMG(pc,b,y);CHKERRQ(ierr);
14138b1f7689SBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGGetNumIterations,(jac->ss_solver,(HYPRE_Int *)&oits));
1414d851a50bSGlenn Hammond   *outits = oits;
1415d851a50bSGlenn Hammond   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
1416d851a50bSGlenn Hammond   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
1417fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,jac->tol));
1418fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,jac->its));
1419d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1420d851a50bSGlenn Hammond }
1421d851a50bSGlenn Hammond 
1422d851a50bSGlenn Hammond 
1423d851a50bSGlenn Hammond #undef __FUNCT__
1424d851a50bSGlenn Hammond #define __FUNCT__ "PCSetUp_SysPFMG"
1425d851a50bSGlenn Hammond PetscErrorCode PCSetUp_SysPFMG(PC pc)
1426d851a50bSGlenn Hammond {
1427d851a50bSGlenn Hammond   PetscErrorCode   ierr;
1428d851a50bSGlenn Hammond   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
1429d851a50bSGlenn Hammond   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct*)(pc->pmat->data);
1430ace3abfcSBarry Smith   PetscBool        flg;
1431d851a50bSGlenn Hammond 
1432d851a50bSGlenn Hammond   PetscFunctionBegin;
1433251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESSTRUCT,&flg);CHKERRQ(ierr);
1434ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESSTRUCT with this preconditioner");
1435d851a50bSGlenn Hammond 
1436d851a50bSGlenn Hammond   /* create the hypre sstruct solver object and set its information */
14372fa5cd67SKarl Rupp   if (ex->ss_solver) PetscStackCallStandard(HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));
1438fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
1439d851a50bSGlenn Hammond   ierr = PCSetFromOptions_SysPFMG(pc);CHKERRQ(ierr);
1440fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetZeroGuess,(ex->ss_solver));
1441fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGSetup,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1442d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1443d851a50bSGlenn Hammond }
1444d851a50bSGlenn Hammond 
1445d851a50bSGlenn Hammond 
1446d851a50bSGlenn Hammond /*MC
1447d851a50bSGlenn Hammond      PCSysPFMG - the hypre SysPFMG multigrid solver
1448d851a50bSGlenn Hammond 
1449d851a50bSGlenn Hammond    Level: advanced
1450d851a50bSGlenn Hammond 
1451d851a50bSGlenn Hammond    Options Database:
1452d851a50bSGlenn Hammond + -pc_syspfmg_its <its> number of iterations of SysPFMG to use as preconditioner
1453d851a50bSGlenn Hammond . -pc_syspfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
1454d851a50bSGlenn Hammond . -pc_syspfmg_num_post_relax <steps> number of smoothing steps after coarse grid
1455d851a50bSGlenn Hammond . -pc_syspfmg_tol <tol> tolerance of SysPFMG
1456d851a50bSGlenn Hammond . -pc_syspfmg_relax_type -relaxation type for the up and down cycles, one of Weighted-Jacobi,Red/Black-Gauss-Seidel
1457d851a50bSGlenn Hammond 
1458d851a50bSGlenn Hammond    Notes:  This is for CELL-centered descretizations
1459d851a50bSGlenn Hammond 
1460f6680f47SSatish Balay            This must be used with the MATHYPRESSTRUCT matrix type.
1461aa219208SBarry Smith            This is less general than in hypre, it supports only one part, and one block per process defined by a PETSc DMDA.
1462d851a50bSGlenn Hammond            Also, only cell-centered variables.
1463d851a50bSGlenn Hammond 
1464d851a50bSGlenn Hammond .seealso:  PCMG, MATHYPRESSTRUCT
1465d851a50bSGlenn Hammond M*/
1466d851a50bSGlenn Hammond 
1467d851a50bSGlenn Hammond #undef __FUNCT__
1468d851a50bSGlenn Hammond #define __FUNCT__ "PCCreate_SysPFMG"
14698cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SysPFMG(PC pc)
1470d851a50bSGlenn Hammond {
1471d851a50bSGlenn Hammond   PetscErrorCode ierr;
1472d851a50bSGlenn Hammond   PC_SysPFMG     *ex;
1473d851a50bSGlenn Hammond 
1474d851a50bSGlenn Hammond   PetscFunctionBegin;
1475*b00a9115SJed Brown   ierr     = PetscNew(&ex);CHKERRQ(ierr); \
1476d851a50bSGlenn Hammond   pc->data = ex;
1477d851a50bSGlenn Hammond 
1478d851a50bSGlenn Hammond   ex->its            = 1;
1479d851a50bSGlenn Hammond   ex->tol            = 1.e-8;
1480d851a50bSGlenn Hammond   ex->relax_type     = 1;
1481d851a50bSGlenn Hammond   ex->num_pre_relax  = 1;
1482d851a50bSGlenn Hammond   ex->num_post_relax = 1;
1483d851a50bSGlenn Hammond 
1484d851a50bSGlenn Hammond   pc->ops->setfromoptions  = PCSetFromOptions_SysPFMG;
1485d851a50bSGlenn Hammond   pc->ops->view            = PCView_SysPFMG;
1486d851a50bSGlenn Hammond   pc->ops->destroy         = PCDestroy_SysPFMG;
1487d851a50bSGlenn Hammond   pc->ops->apply           = PCApply_SysPFMG;
1488d851a50bSGlenn Hammond   pc->ops->applyrichardson = PCApplyRichardson_SysPFMG;
1489d851a50bSGlenn Hammond   pc->ops->setup           = PCSetUp_SysPFMG;
14902fa5cd67SKarl Rupp 
1491ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ex->hcomm));CHKERRQ(ierr);
1492fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
1493d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1494d851a50bSGlenn Hammond }
1495