xref: /petsc/src/ksp/pc/impls/hypre/hypre.c (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
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 
8*b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>          /*I "petscpc.h" I*/
9c6db04a5SJed Brown #include <../src/dm/impls/da/hypre/mhyp.h>
1016d9e3a6SLisandro Dalcin 
1116d9e3a6SLisandro Dalcin /*
1216d9e3a6SLisandro Dalcin    Private context (data structure) for the  preconditioner.
1316d9e3a6SLisandro Dalcin */
1416d9e3a6SLisandro Dalcin typedef struct {
1516d9e3a6SLisandro Dalcin   HYPRE_Solver       hsolver;
1616d9e3a6SLisandro Dalcin   HYPRE_IJMatrix     ij;
1716d9e3a6SLisandro Dalcin   HYPRE_IJVector     b,x;
1816d9e3a6SLisandro Dalcin 
194ddd07fcSJed Brown   HYPRE_Int          (*destroy)(HYPRE_Solver);
204ddd07fcSJed Brown   HYPRE_Int          (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
214ddd07fcSJed Brown   HYPRE_Int          (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
2216d9e3a6SLisandro Dalcin 
2316d9e3a6SLisandro Dalcin   MPI_Comm           comm_hypre;
2416d9e3a6SLisandro Dalcin   char              *hypre_type;
2516d9e3a6SLisandro Dalcin 
2616d9e3a6SLisandro Dalcin   /* options for Pilut and BoomerAMG*/
274ddd07fcSJed Brown   PetscInt           maxiter;
2816d9e3a6SLisandro Dalcin   double             tol;
2916d9e3a6SLisandro Dalcin 
3016d9e3a6SLisandro Dalcin   /* options for Pilut */
314ddd07fcSJed Brown   PetscInt           factorrowsize;
3216d9e3a6SLisandro Dalcin 
3316d9e3a6SLisandro Dalcin   /* options for ParaSails */
344ddd07fcSJed Brown   PetscInt           nlevels;
3516d9e3a6SLisandro Dalcin   double             threshhold;
3616d9e3a6SLisandro Dalcin   double             filter;
374ddd07fcSJed Brown   PetscInt           sym;
3816d9e3a6SLisandro Dalcin   double             loadbal;
394ddd07fcSJed Brown   PetscInt           logging;
404ddd07fcSJed Brown   PetscInt           ruse;
414ddd07fcSJed Brown   PetscInt           symt;
4216d9e3a6SLisandro Dalcin 
4316d9e3a6SLisandro Dalcin   /* options for Euclid */
44ace3abfcSBarry Smith   PetscBool          bjilu;
454ddd07fcSJed Brown   PetscInt           levels;
4616d9e3a6SLisandro Dalcin 
4716d9e3a6SLisandro Dalcin   /* options for Euclid and BoomerAMG */
48ace3abfcSBarry Smith   PetscBool          printstatistics;
4916d9e3a6SLisandro Dalcin 
5016d9e3a6SLisandro Dalcin   /* options for BoomerAMG */
514ddd07fcSJed Brown   PetscInt           cycletype;
524ddd07fcSJed Brown   PetscInt           maxlevels;
5316d9e3a6SLisandro Dalcin   double             strongthreshold;
5416d9e3a6SLisandro Dalcin   double             maxrowsum;
554ddd07fcSJed Brown   PetscInt           gridsweeps[3];
564ddd07fcSJed Brown   PetscInt           coarsentype;
574ddd07fcSJed Brown   PetscInt           measuretype;
584ddd07fcSJed Brown   PetscInt           relaxtype[3];
5916d9e3a6SLisandro Dalcin   double             relaxweight;
6016d9e3a6SLisandro Dalcin   double             outerrelaxweight;
614ddd07fcSJed Brown   PetscInt           relaxorder;
6216d9e3a6SLisandro Dalcin   double             truncfactor;
63ace3abfcSBarry Smith   PetscBool          applyrichardson;
644ddd07fcSJed Brown   PetscInt           pmax;
654ddd07fcSJed Brown   PetscInt           interptype;
664ddd07fcSJed Brown   PetscInt           agg_nl;
674ddd07fcSJed Brown   PetscInt           agg_num_paths;
684ddd07fcSJed Brown   PetscInt           nodal_coarsen;
69ace3abfcSBarry Smith   PetscBool          nodal_relax;
704ddd07fcSJed Brown   PetscInt           nodal_relax_levels;
7116d9e3a6SLisandro Dalcin } PC_HYPRE;
7216d9e3a6SLisandro Dalcin 
7316d9e3a6SLisandro Dalcin 
7416d9e3a6SLisandro Dalcin #undef __FUNCT__
7516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetUp_HYPRE"
7616d9e3a6SLisandro Dalcin static PetscErrorCode PCSetUp_HYPRE(PC pc)
7716d9e3a6SLisandro Dalcin {
7816d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
7916d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
8016d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
8116d9e3a6SLisandro Dalcin   HYPRE_ParVector    bv,xv;
8216d9e3a6SLisandro Dalcin   PetscInt           bs;
8316d9e3a6SLisandro Dalcin 
8416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
8516d9e3a6SLisandro Dalcin   if (!jac->hypre_type) {
8602a17cd4SBarry Smith     ierr = PCHYPRESetType(pc,"boomeramg");CHKERRQ(ierr);
8716d9e3a6SLisandro Dalcin   }
885f5c5b43SBarry Smith 
895f5c5b43SBarry Smith   if (pc->setupcalled) {
905f5c5b43SBarry Smith     /* always destroy the old matrix and create a new memory;
915f5c5b43SBarry Smith        hope this does not churn the memory too much. The problem
925f5c5b43SBarry Smith        is I do not know if it is possible to put the matrix back to
935f5c5b43SBarry Smith        its initial state so that we can directly copy the values
945f5c5b43SBarry Smith        the second time through. */
9530e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_IJMatrixDestroy,(jac->ij));
965f5c5b43SBarry Smith     jac->ij = 0;
9716d9e3a6SLisandro Dalcin   }
985f5c5b43SBarry Smith 
9916d9e3a6SLisandro Dalcin   if (!jac->ij) { /* create the matrix the first time through */
10016d9e3a6SLisandro Dalcin     ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr);
10116d9e3a6SLisandro Dalcin   }
10216d9e3a6SLisandro Dalcin   if (!jac->b) { /* create the vectors the first time through */
10316d9e3a6SLisandro Dalcin     Vec x,b;
10416d9e3a6SLisandro Dalcin     ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr);
10516d9e3a6SLisandro Dalcin     ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr);
10616d9e3a6SLisandro Dalcin     ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr);
1076bf464f9SBarry Smith     ierr = VecDestroy(&x);CHKERRQ(ierr);
1086bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
10916d9e3a6SLisandro Dalcin   }
1105f5c5b43SBarry Smith 
11116d9e3a6SLisandro Dalcin   /* special case for BoomerAMG */
11216d9e3a6SLisandro Dalcin   if (jac->setup == HYPRE_BoomerAMGSetup) {
11316d9e3a6SLisandro Dalcin     ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr);
11416d9e3a6SLisandro Dalcin     if (bs > 1) {
11530e6f737SJed Brown       PetscStackCallHypre(0,HYPRE_BoomerAMGSetNumFunctions,(jac->hsolver,bs));
11616d9e3a6SLisandro Dalcin     }
11716d9e3a6SLisandro Dalcin   };
11816d9e3a6SLisandro Dalcin   ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr);
11930e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
12030e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJVectorGetObject,(jac->b,(void**)&bv));
12130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJVectorGetObject,(jac->x,(void**)&xv));
12230e6f737SJed Brown   PetscStackCallHypre("HYPRE_SetupXXX",(*jac->setup),(jac->hsolver,hmat,bv,xv));
12316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
12416d9e3a6SLisandro Dalcin }
12516d9e3a6SLisandro Dalcin 
12616d9e3a6SLisandro Dalcin /*
12716d9e3a6SLisandro Dalcin     Replaces the address where the HYPRE vector points to its data with the address of
12816d9e3a6SLisandro Dalcin   PETSc's data. Saves the old address so it can be reset when we are finished with it.
12916d9e3a6SLisandro Dalcin   Allows use to get the data into a HYPRE vector without the cost of memcopies
13016d9e3a6SLisandro Dalcin */
13116d9e3a6SLisandro Dalcin #define HYPREReplacePointer(b,newvalue,savedvalue) {\
13216d9e3a6SLisandro Dalcin    hypre_ParVector *par_vector   = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));\
13316d9e3a6SLisandro Dalcin    hypre_Vector    *local_vector = hypre_ParVectorLocalVector(par_vector);\
13416d9e3a6SLisandro Dalcin    savedvalue         = local_vector->data;\
13516d9e3a6SLisandro Dalcin    local_vector->data = newvalue;}
13616d9e3a6SLisandro Dalcin 
13716d9e3a6SLisandro Dalcin #undef __FUNCT__
13816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApply_HYPRE"
13916d9e3a6SLisandro Dalcin static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
14016d9e3a6SLisandro Dalcin {
14116d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
14216d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
14316d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
14416d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
14516d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
14616d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
1474ddd07fcSJed Brown   PetscInt           hierr;
14816d9e3a6SLisandro Dalcin 
14916d9e3a6SLisandro Dalcin   PetscFunctionBegin;
15016d9e3a6SLisandro Dalcin   if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);}
15116d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
15216d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
15316d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
15416d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
15516d9e3a6SLisandro Dalcin 
15630e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
15730e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJVectorGetObject,(jac->b,(void**)&jbv));
15830e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJVectorGetObject,(jac->x,(void**)&jxv));
15916d9e3a6SLisandro Dalcin   hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv);
16065e19b50SBarry Smith   if (hierr && hierr != HYPRE_ERROR_CONV) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
16116d9e3a6SLisandro Dalcin   if (hierr) hypre__global_error = 0;
16216d9e3a6SLisandro Dalcin 
16316d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
16416d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
16516d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
16616d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
16716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
16816d9e3a6SLisandro Dalcin }
16916d9e3a6SLisandro Dalcin 
17016d9e3a6SLisandro Dalcin #undef __FUNCT__
17116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCDestroy_HYPRE"
17216d9e3a6SLisandro Dalcin static PetscErrorCode PCDestroy_HYPRE(PC pc)
17316d9e3a6SLisandro Dalcin {
17416d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
17516d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
17616d9e3a6SLisandro Dalcin 
17716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
17830e6f737SJed Brown   if (jac->ij) PetscStackCallHypre(0,HYPRE_IJMatrixDestroy,(jac->ij));
17930e6f737SJed Brown   if (jac->b)  PetscStackCallHypre(0,HYPRE_IJVectorDestroy,(jac->b));
18030e6f737SJed Brown   if (jac->x)  PetscStackCallHypre(0,HYPRE_IJVectorDestroy,(jac->x));
18130e6f737SJed Brown   if (jac->destroy) PetscStackCallHypre("HYPRE_DistroyXXX",(*jac->destroy),(jac->hsolver));
182503cfb0cSBarry Smith   ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr);
18316d9e3a6SLisandro Dalcin   if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);}
184c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
18516d9e3a6SLisandro Dalcin 
18616d9e3a6SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr);
18716d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);CHKERRQ(ierr);
18816d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);CHKERRQ(ierr);
18916d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
19016d9e3a6SLisandro Dalcin }
19116d9e3a6SLisandro Dalcin 
19216d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
19316d9e3a6SLisandro Dalcin #undef __FUNCT__
19416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut"
19516d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc)
19616d9e3a6SLisandro Dalcin {
19716d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
19816d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
199ace3abfcSBarry Smith   PetscBool      flag;
20016d9e3a6SLisandro Dalcin 
20116d9e3a6SLisandro Dalcin   PetscFunctionBegin;
20216d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr);
20316d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr);
20430e6f737SJed Brown   if (flag) PetscStackCallHypre(0,HYPRE_ParCSRPilutSetMaxIter,(jac->hsolver,jac->maxiter));
20516d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr);
20630e6f737SJed Brown   if (flag) PetscStackCallHypre(0,HYPRE_ParCSRPilutSetDropTolerance,(jac->hsolver,jac->tol));
20716d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr);
20830e6f737SJed Brown   if (flag) PetscStackCallHypre(0,HYPRE_ParCSRPilutSetFactorRowSize,(jac->hsolver,jac->factorrowsize));
20916d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
21016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
21116d9e3a6SLisandro Dalcin }
21216d9e3a6SLisandro Dalcin 
21316d9e3a6SLisandro Dalcin #undef __FUNCT__
21416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Pilut"
21516d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
21616d9e3a6SLisandro Dalcin {
21716d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
21816d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
219ace3abfcSBarry Smith   PetscBool      iascii;
22016d9e3a6SLisandro Dalcin 
22116d9e3a6SLisandro Dalcin   PetscFunctionBegin;
2222692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
22316d9e3a6SLisandro Dalcin   if (iascii) {
22416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");CHKERRQ(ierr);
22516d9e3a6SLisandro Dalcin     if (jac->maxiter != PETSC_DEFAULT) {
22616d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr);
22716d9e3a6SLisandro Dalcin     } else {
22816d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr);
22916d9e3a6SLisandro Dalcin     }
23016d9e3a6SLisandro Dalcin     if (jac->tol != PETSC_DEFAULT) {
23116d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr);
23216d9e3a6SLisandro Dalcin     } else {
23316d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr);
23416d9e3a6SLisandro Dalcin     }
23516d9e3a6SLisandro Dalcin     if (jac->factorrowsize != PETSC_DEFAULT) {
23616d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr);
23716d9e3a6SLisandro Dalcin     } else {
23816d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");CHKERRQ(ierr);
23916d9e3a6SLisandro Dalcin     }
24016d9e3a6SLisandro Dalcin   }
24116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
24216d9e3a6SLisandro Dalcin }
24316d9e3a6SLisandro Dalcin 
24416d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
24516d9e3a6SLisandro Dalcin #undef __FUNCT__
24616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid"
24716d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc)
24816d9e3a6SLisandro Dalcin {
24916d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
25016d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
251ace3abfcSBarry Smith   PetscBool      flag;
252390e7148SBarry Smith   char           *args[8],levels[16];
253390e7148SBarry Smith   PetscInt       cnt = 0;
25416d9e3a6SLisandro Dalcin 
25516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
25616d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr);
25716d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr);
25816d9e3a6SLisandro Dalcin   if (flag) {
25965e19b50SBarry Smith     if (jac->levels < 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels);
2604ddd07fcSJed Brown     ierr = PetscSNPrintf(levels,sizeof levels,"%D",jac->levels);CHKERRQ(ierr);
261390e7148SBarry Smith     args[cnt++] = (char*)"-level"; args[cnt++] = levels;
26216d9e3a6SLisandro Dalcin   }
263acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);CHKERRQ(ierr);
26416d9e3a6SLisandro Dalcin   if (jac->bjilu) {
265390e7148SBarry Smith     args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1";
26616d9e3a6SLisandro Dalcin   }
26716d9e3a6SLisandro Dalcin 
268acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);CHKERRQ(ierr);
26916d9e3a6SLisandro Dalcin   if (jac->printstatistics) {
270390e7148SBarry Smith     args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1";
271390e7148SBarry Smith     args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1";
27216d9e3a6SLisandro Dalcin   }
27316d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
27430e6f737SJed Brown   if (cnt) PetscStackCallHypre(0,HYPRE_EuclidSetParams,(jac->hsolver,cnt,args));
27516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
27616d9e3a6SLisandro Dalcin }
27716d9e3a6SLisandro Dalcin 
27816d9e3a6SLisandro Dalcin #undef __FUNCT__
27916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Euclid"
28016d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
28116d9e3a6SLisandro Dalcin {
28216d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
28316d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
284ace3abfcSBarry Smith   PetscBool      iascii;
28516d9e3a6SLisandro Dalcin 
28616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
2872692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
28816d9e3a6SLisandro Dalcin   if (iascii) {
28916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n");CHKERRQ(ierr);
29016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr);
29116d9e3a6SLisandro Dalcin     if (jac->bjilu) {
29216d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr);
29316d9e3a6SLisandro Dalcin     }
29416d9e3a6SLisandro Dalcin   }
29516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
29616d9e3a6SLisandro Dalcin }
29716d9e3a6SLisandro Dalcin 
29816d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
29916d9e3a6SLisandro Dalcin 
30016d9e3a6SLisandro Dalcin #undef __FUNCT__
30116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG"
30216d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
30316d9e3a6SLisandro Dalcin {
30416d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
30516d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
30616d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
30716d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
30816d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
30916d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
3104ddd07fcSJed Brown   PetscInt           hierr;
31116d9e3a6SLisandro Dalcin 
31216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
31316d9e3a6SLisandro Dalcin   ierr = VecSet(x,0.0);CHKERRQ(ierr);
31416d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
31516d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
31616d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
31716d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
31816d9e3a6SLisandro Dalcin 
31930e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat));
32030e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJVectorGetObject,(jac->b,(void**)&jbv));
32130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_IJVectorGetObject,(jac->x,(void**)&jxv));
32216d9e3a6SLisandro Dalcin 
32316d9e3a6SLisandro Dalcin   hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv);
32416d9e3a6SLisandro Dalcin   /* error code of 1 in BoomerAMG merely means convergence not achieved */
325e32f2f54SBarry Smith   if (hierr && (hierr != 1)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
32616d9e3a6SLisandro Dalcin   if (hierr) hypre__global_error = 0;
32716d9e3a6SLisandro Dalcin 
32816d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
32916d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
33016d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
33116d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
33216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
33316d9e3a6SLisandro Dalcin }
33416d9e3a6SLisandro Dalcin 
33516d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
3360f1074feSSatish Balay static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"};
33716d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
3380f1074feSSatish Balay static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","","SOR/Jacobi","backward-SOR/Jacobi","","symmetric-SOR/Jacobi",
3390f1074feSSatish Balay                                                   "","","Gaussian-elimination"};
3400f1074feSSatish Balay static const char *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i",
3410f1074feSSatish Balay                                                   "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"};
34216d9e3a6SLisandro Dalcin #undef __FUNCT__
34316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG"
34416d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
34516d9e3a6SLisandro Dalcin {
34616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
34716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
3484ddd07fcSJed Brown   PetscInt       n,indx,level;
349ace3abfcSBarry Smith   PetscBool      flg, tmp_truth;
35016d9e3a6SLisandro Dalcin   double         tmpdbl, twodbl[2];
35116d9e3a6SLisandro Dalcin 
35216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
35316d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr);
3544336a9eeSBarry Smith   ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType+1,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr);
35516d9e3a6SLisandro Dalcin   if (flg) {
3564336a9eeSBarry Smith     jac->cycletype = indx+1;
35730e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
35816d9e3a6SLisandro Dalcin   }
35916d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr);
36016d9e3a6SLisandro Dalcin   if (flg) {
36165e19b50SBarry Smith     if (jac->maxlevels < 2) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
36230e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
36316d9e3a6SLisandro Dalcin   }
36416d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr);
36516d9e3a6SLisandro Dalcin   if (flg) {
36665e19b50SBarry Smith     if (jac->maxiter < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
36730e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
36816d9e3a6SLisandro Dalcin   }
3690f1074feSSatish 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);
37016d9e3a6SLisandro Dalcin   if (flg) {
37165e19b50SBarry Smith     if (jac->tol < 0.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be greater than or equal to zero",jac->tol);
37230e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
37316d9e3a6SLisandro Dalcin   }
37416d9e3a6SLisandro Dalcin 
3750f1074feSSatish Balay   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr);
37616d9e3a6SLisandro Dalcin   if (flg) {
37765e19b50SBarry Smith     if (jac->truncfactor < 0.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor);
37830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
37916d9e3a6SLisandro Dalcin   }
38016d9e3a6SLisandro Dalcin 
3810f1074feSSatish 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);
3820f1074feSSatish Balay   if (flg) {
38365e19b50SBarry Smith     if (jac->pmax < 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"P_max %G must be greater than or equal to zero",jac->pmax);
38430e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
3850f1074feSSatish Balay   }
3860f1074feSSatish Balay 
3870f1074feSSatish Balay  ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);CHKERRQ(ierr);
3880f1074feSSatish Balay   if (flg) {
38965e19b50SBarry Smith      if (jac->agg_nl < 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %G must be greater than or equal to zero",jac->agg_nl);
3900f1074feSSatish Balay 
39130e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
3920f1074feSSatish Balay   }
3930f1074feSSatish Balay 
3940f1074feSSatish Balay 
3950f1074feSSatish 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);
3960f1074feSSatish Balay   if (flg) {
39765e19b50SBarry Smith      if (jac->agg_num_paths < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of paths %G must be greater than or equal to 1",jac->agg_num_paths);
3980f1074feSSatish Balay 
39930e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
4000f1074feSSatish Balay   }
4010f1074feSSatish Balay 
4020f1074feSSatish Balay 
40316d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr);
40416d9e3a6SLisandro Dalcin   if (flg) {
40565e19b50SBarry Smith     if (jac->strongthreshold < 0.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold);
40630e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
40716d9e3a6SLisandro Dalcin   }
40816d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr);
40916d9e3a6SLisandro Dalcin   if (flg) {
41065e19b50SBarry Smith     if (jac->maxrowsum < 0.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum);
41165e19b50SBarry Smith     if (jac->maxrowsum > 1.0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum);
41230e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
41316d9e3a6SLisandro Dalcin   }
41416d9e3a6SLisandro Dalcin 
41516d9e3a6SLisandro Dalcin   /* Grid sweeps */
4160f1074feSSatish 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);
41716d9e3a6SLisandro Dalcin   if (flg) {
41830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver,indx));
41916d9e3a6SLisandro Dalcin     /* modify the jac structure so we can view the updated options with PC_View */
42016d9e3a6SLisandro Dalcin     jac->gridsweeps[0] = indx;
4210f1074feSSatish Balay     jac->gridsweeps[1] = indx;
4220f1074feSSatish Balay     /*defaults coarse to 1 */
4230f1074feSSatish Balay     jac->gridsweeps[2] = 1;
42416d9e3a6SLisandro Dalcin   }
4250f1074feSSatish Balay 
4260f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx ,&flg);CHKERRQ(ierr);
42716d9e3a6SLisandro Dalcin   if (flg) {
42830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 1));
4290f1074feSSatish Balay     jac->gridsweeps[0] = indx;
43016d9e3a6SLisandro Dalcin   }
43116d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr);
43216d9e3a6SLisandro Dalcin   if (flg) {
43330e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 2));
4340f1074feSSatish Balay     jac->gridsweeps[1] = indx;
43516d9e3a6SLisandro Dalcin   }
4360f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr);
43716d9e3a6SLisandro Dalcin   if (flg) {
43830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 3));
4390f1074feSSatish Balay     jac->gridsweeps[2] = indx;
44016d9e3a6SLisandro Dalcin   }
44116d9e3a6SLisandro Dalcin 
44216d9e3a6SLisandro Dalcin   /* Relax type */
4430f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_all","Relax type for the up and down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
44416d9e3a6SLisandro Dalcin   if (flg) {
4450f1074feSSatish Balay     jac->relaxtype[0] = jac->relaxtype[1]  = indx;
44630e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, indx));
4470f1074feSSatish Balay     /* by default, coarse type set to 9 */
4480f1074feSSatish Balay     jac->relaxtype[2] = 9;
4490f1074feSSatish Balay 
45016d9e3a6SLisandro Dalcin   }
4510f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
45216d9e3a6SLisandro Dalcin   if (flg) {
45316d9e3a6SLisandro Dalcin     jac->relaxtype[0] = indx;
45430e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 1));
45516d9e3a6SLisandro Dalcin   }
4560f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
45716d9e3a6SLisandro Dalcin   if (flg) {
4580f1074feSSatish Balay     jac->relaxtype[1] = indx;
45930e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 2));
46016d9e3a6SLisandro Dalcin   }
46116d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr);
46216d9e3a6SLisandro Dalcin   if (flg) {
4630f1074feSSatish Balay     jac->relaxtype[2] = indx;
46430e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 3));
46516d9e3a6SLisandro Dalcin   }
46616d9e3a6SLisandro Dalcin 
46716d9e3a6SLisandro Dalcin   /* Relaxation Weight */
46816d9e3a6SLisandro 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);
46916d9e3a6SLisandro Dalcin   if (flg) {
47030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetRelaxWt,(jac->hsolver,tmpdbl));
47116d9e3a6SLisandro Dalcin     jac->relaxweight = tmpdbl;
47216d9e3a6SLisandro Dalcin   }
47316d9e3a6SLisandro Dalcin 
47416d9e3a6SLisandro Dalcin   n=2;
47516d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
47616d9e3a6SLisandro 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);
47716d9e3a6SLisandro Dalcin   if (flg) {
47816d9e3a6SLisandro Dalcin     if (n == 2) {
47916d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
48030e6f737SJed Brown       PetscStackCallHypre(0,HYPRE_BoomerAMGSetLevelRelaxWt,(jac->hsolver,twodbl[0],indx));
48116d9e3a6SLisandro Dalcin     } else {
48265e19b50SBarry Smith       SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relax weight level: you must provide 2 values separated by a comma (and no space), you provided %d",n);
48316d9e3a6SLisandro Dalcin     }
48416d9e3a6SLisandro Dalcin   }
48516d9e3a6SLisandro Dalcin 
48616d9e3a6SLisandro Dalcin   /* Outer relaxation Weight */
48716d9e3a6SLisandro 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);
48816d9e3a6SLisandro Dalcin   if (flg) {
48930e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetOuterWt,( jac->hsolver, tmpdbl));
49016d9e3a6SLisandro Dalcin     jac->outerrelaxweight = tmpdbl;
49116d9e3a6SLisandro Dalcin   }
49216d9e3a6SLisandro Dalcin 
49316d9e3a6SLisandro Dalcin   n=2;
49416d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
49516d9e3a6SLisandro 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);
49616d9e3a6SLisandro Dalcin   if (flg) {
49716d9e3a6SLisandro Dalcin     if (n == 2) {
49816d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
49930e6f737SJed Brown       PetscStackCallHypre(0,HYPRE_BoomerAMGSetLevelOuterWt,( jac->hsolver, twodbl[0], indx));
50016d9e3a6SLisandro Dalcin     } else {
50165e19b50SBarry Smith       SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relax weight outer level: You must provide 2 values separated by a comma (and no space), you provided %d",n);
50216d9e3a6SLisandro Dalcin     }
50316d9e3a6SLisandro Dalcin   }
50416d9e3a6SLisandro Dalcin 
50516d9e3a6SLisandro Dalcin   /* the Relax Order */
506acfcf0e5SJed Brown   ierr = PetscOptionsBool( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
50716d9e3a6SLisandro Dalcin 
50816d9e3a6SLisandro Dalcin   if (flg) {
50916d9e3a6SLisandro Dalcin     jac->relaxorder = 0;
51030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
51116d9e3a6SLisandro Dalcin   }
51216d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr);
51316d9e3a6SLisandro Dalcin   if (flg) {
51416d9e3a6SLisandro Dalcin     jac->measuretype = indx;
51530e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
51616d9e3a6SLisandro Dalcin   }
5170f1074feSSatish Balay   /* update list length 3/07 */
5180f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,11,HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr);
51916d9e3a6SLisandro Dalcin   if (flg) {
52016d9e3a6SLisandro Dalcin     jac->coarsentype = indx;
52130e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
52216d9e3a6SLisandro Dalcin   }
5230f1074feSSatish Balay 
5240f1074feSSatish Balay   /* new 3/07 */
5250f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,14,HYPREBoomerAMGInterpType[0],&indx,&flg);CHKERRQ(ierr);
5260f1074feSSatish Balay   if (flg) {
5270f1074feSSatish Balay     jac->interptype = indx;
52830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
5290f1074feSSatish Balay   }
5300f1074feSSatish Balay 
531b96a4a96SBarry Smith   ierr = PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);CHKERRQ(ierr);
53216d9e3a6SLisandro Dalcin   if (flg) {
533b96a4a96SBarry Smith     level = 3;
53416d9e3a6SLisandro Dalcin     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
535b96a4a96SBarry Smith     jac->printstatistics = PETSC_TRUE;
53630e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetPrintLevel,(jac->hsolver,level));
5372ae77aedSBarry Smith   }
5382ae77aedSBarry Smith 
539b96a4a96SBarry Smith   ierr = PetscOptionsName("-pc_hypre_boomeramg_print_debug","Print debug information","None",&flg);CHKERRQ(ierr);
5402ae77aedSBarry Smith   if (flg) {
541b96a4a96SBarry Smith     level = 3;
5422ae77aedSBarry Smith     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
543b96a4a96SBarry Smith     jac->printstatistics = PETSC_TRUE;
54430e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetDebugFlag,(jac->hsolver,level));
54516d9e3a6SLisandro Dalcin   }
5468f87f92bSBarry Smith 
547acfcf0e5SJed Brown   ierr = PetscOptionsBool( "-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5488f87f92bSBarry Smith   if (flg && tmp_truth) {
5498f87f92bSBarry Smith     jac->nodal_coarsen = 1;
55030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetNodal,(jac->hsolver,1));
5518f87f92bSBarry Smith   }
5528f87f92bSBarry Smith 
553acfcf0e5SJed Brown   ierr = PetscOptionsBool( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5548f87f92bSBarry Smith   if (flg && tmp_truth) {
5558f87f92bSBarry Smith     PetscInt tmp_int;
5568f87f92bSBarry Smith     ierr = PetscOptionsInt( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);CHKERRQ(ierr);
5578f87f92bSBarry Smith     if (flg) jac->nodal_relax_levels = tmp_int;
55830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetSmoothType,(jac->hsolver,6));
55930e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetDomainType,(jac->hsolver,1));
56030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetOverlap,(jac->hsolver,0));
56130e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetSmoothNumLevels,(jac->hsolver,jac->nodal_relax_levels));
5628f87f92bSBarry Smith   }
5638f87f92bSBarry Smith 
56416d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
56516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
56616d9e3a6SLisandro Dalcin }
56716d9e3a6SLisandro Dalcin 
56816d9e3a6SLisandro Dalcin #undef __FUNCT__
56916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG"
570ace3abfcSBarry 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)
57116d9e3a6SLisandro Dalcin {
57216d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
57316d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
5744ddd07fcSJed Brown   PetscInt       oits;
57516d9e3a6SLisandro Dalcin 
57616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
57730e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,its*jac->maxiter));
57830e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_BoomerAMGSetTol,(jac->hsolver,rtol));
57916d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_TRUE;
58016d9e3a6SLisandro Dalcin   ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr);
58116d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_FALSE;
58230e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_BoomerAMGGetNumIterations,(jac->hsolver,&oits));
5834d0a8057SBarry Smith   *outits = oits;
5844d0a8057SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
5854d0a8057SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
58630e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
58730e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
58816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
58916d9e3a6SLisandro Dalcin }
59016d9e3a6SLisandro Dalcin 
59116d9e3a6SLisandro Dalcin 
59216d9e3a6SLisandro Dalcin #undef __FUNCT__
59316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_BoomerAMG"
59416d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
59516d9e3a6SLisandro Dalcin {
59616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
59716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
598ace3abfcSBarry Smith   PetscBool      iascii;
59916d9e3a6SLisandro Dalcin 
60016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
6012692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
60216d9e3a6SLisandro Dalcin   if (iascii) {
60316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr);
60416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr);
60516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr);
60616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr);
60716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr);
60816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr);
6090f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation truncation factor %G\n",jac->truncfactor);CHKERRQ(ierr);
6100f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);CHKERRQ(ierr);
6110f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);CHKERRQ(ierr);
6120f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);CHKERRQ(ierr);
6130f1074feSSatish Balay 
61416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr);
61516d9e3a6SLisandro Dalcin 
6160f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps down         %d\n",jac->gridsweeps[0]);CHKERRQ(ierr);
6170f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps up           %d\n",jac->gridsweeps[1]);CHKERRQ(ierr);
6180f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on coarse    %d\n",jac->gridsweeps[2]);CHKERRQ(ierr);
61916d9e3a6SLisandro Dalcin 
6200f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr);
6210f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr);
6220f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr);
62316d9e3a6SLisandro Dalcin 
62416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax weight  (all)      %G\n",jac->relaxweight);CHKERRQ(ierr);
62516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr);
62616d9e3a6SLisandro Dalcin 
62716d9e3a6SLisandro Dalcin     if (jac->relaxorder) {
62816d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr);
62916d9e3a6SLisandro Dalcin     } else {
63016d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr);
63116d9e3a6SLisandro Dalcin     }
63216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr);
63316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr);
6340f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation type  %s\n",HYPREBoomerAMGInterpType[jac->interptype]);CHKERRQ(ierr);
6358f87f92bSBarry Smith     if (jac->nodal_coarsen) {
6368f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal())\n");CHKERRQ(ierr);
6378f87f92bSBarry Smith     }
6388f87f92bSBarry Smith     if (jac->nodal_relax) {
6398f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);CHKERRQ(ierr);
6408f87f92bSBarry Smith     }
64116d9e3a6SLisandro Dalcin   }
64216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
64316d9e3a6SLisandro Dalcin }
64416d9e3a6SLisandro Dalcin 
64516d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
64616d9e3a6SLisandro Dalcin #undef __FUNCT__
64716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails"
64816d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc)
64916d9e3a6SLisandro Dalcin {
65016d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
65116d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
6524ddd07fcSJed Brown   PetscInt       indx;
653ace3abfcSBarry Smith   PetscBool      flag;
65416d9e3a6SLisandro Dalcin   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};
65516d9e3a6SLisandro Dalcin 
65616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
65716d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr);
65816d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr);
65916d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr);
66016d9e3a6SLisandro Dalcin   if (flag) {
66130e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels));
66216d9e3a6SLisandro Dalcin   }
66316d9e3a6SLisandro Dalcin 
66416d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr);
66516d9e3a6SLisandro Dalcin   if (flag) {
66630e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter));
66716d9e3a6SLisandro Dalcin   }
66816d9e3a6SLisandro Dalcin 
66916d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr);
67016d9e3a6SLisandro Dalcin   if (flag) {
67130e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal));
67216d9e3a6SLisandro Dalcin   }
67316d9e3a6SLisandro Dalcin 
674acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_parasails_logging","Print info to screen","None",(PetscBool)jac->logging,(PetscBool*)&jac->logging,&flag);CHKERRQ(ierr);
67516d9e3a6SLisandro Dalcin   if (flag) {
67630e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging));
67716d9e3a6SLisandro Dalcin   }
67816d9e3a6SLisandro Dalcin 
679acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscBool)jac->ruse,(PetscBool*)&jac->ruse,&flag);CHKERRQ(ierr);
68016d9e3a6SLisandro Dalcin   if (flag) {
68130e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse));
68216d9e3a6SLisandro Dalcin   }
68316d9e3a6SLisandro Dalcin 
68416d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);CHKERRQ(ierr);
68516d9e3a6SLisandro Dalcin   if (flag) {
68616d9e3a6SLisandro Dalcin     jac->symt = indx;
68730e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt));
68816d9e3a6SLisandro Dalcin   }
68916d9e3a6SLisandro Dalcin 
69016d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
69116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
69216d9e3a6SLisandro Dalcin }
69316d9e3a6SLisandro Dalcin 
69416d9e3a6SLisandro Dalcin #undef __FUNCT__
69516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_ParaSails"
69616d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
69716d9e3a6SLisandro Dalcin {
69816d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
69916d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
700ace3abfcSBarry Smith   PetscBool      iascii;
70116d9e3a6SLisandro Dalcin   const char     *symt = 0;;
70216d9e3a6SLisandro Dalcin 
70316d9e3a6SLisandro Dalcin   PetscFunctionBegin;
7042692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
70516d9e3a6SLisandro Dalcin   if (iascii) {
70616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");CHKERRQ(ierr);
70716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr);
70816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr);
70916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr);
71016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr);
711ace3abfcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscBools[jac->ruse]);CHKERRQ(ierr);
712ace3abfcSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscBools[jac->logging]);CHKERRQ(ierr);
71316d9e3a6SLisandro Dalcin     if (!jac->symt) {
71416d9e3a6SLisandro Dalcin       symt = "nonsymmetric matrix and preconditioner";
71516d9e3a6SLisandro Dalcin     } else if (jac->symt == 1) {
71616d9e3a6SLisandro Dalcin       symt = "SPD matrix and preconditioner";
71716d9e3a6SLisandro Dalcin     } else if (jac->symt == 2) {
71816d9e3a6SLisandro Dalcin       symt = "nonsymmetric matrix but SPD preconditioner";
71916d9e3a6SLisandro Dalcin     } else {
72065e19b50SBarry Smith       SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
72116d9e3a6SLisandro Dalcin     }
72216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr);
72316d9e3a6SLisandro Dalcin   }
72416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
72516d9e3a6SLisandro Dalcin }
72616d9e3a6SLisandro Dalcin /* ---------------------------------------------------------------------------------*/
72716d9e3a6SLisandro Dalcin 
72816d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
72916d9e3a6SLisandro Dalcin #undef __FUNCT__
73016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType_HYPRE"
7317087cfbeSBarry Smith PetscErrorCode  PCHYPREGetType_HYPRE(PC pc,const char *name[])
73216d9e3a6SLisandro Dalcin {
73316d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
73416d9e3a6SLisandro Dalcin 
73516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
73616d9e3a6SLisandro Dalcin   *name = jac->hypre_type;
73716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
73816d9e3a6SLisandro Dalcin }
73916d9e3a6SLisandro Dalcin EXTERN_C_END
74016d9e3a6SLisandro Dalcin 
74116d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
74216d9e3a6SLisandro Dalcin #undef __FUNCT__
74316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType_HYPRE"
7447087cfbeSBarry Smith PetscErrorCode  PCHYPRESetType_HYPRE(PC pc,const char name[])
74516d9e3a6SLisandro Dalcin {
74616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
74716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
748ace3abfcSBarry Smith   PetscBool      flag;
74916d9e3a6SLisandro Dalcin 
75016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
75116d9e3a6SLisandro Dalcin   if (jac->hypre_type) {
75216d9e3a6SLisandro Dalcin     ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr);
753e7e72b3dSBarry Smith     if (!flag) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
75416d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
75516d9e3a6SLisandro Dalcin   } else {
75616d9e3a6SLisandro Dalcin     ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr);
75716d9e3a6SLisandro Dalcin   }
75816d9e3a6SLisandro Dalcin 
75916d9e3a6SLisandro Dalcin   jac->maxiter            = PETSC_DEFAULT;
76016d9e3a6SLisandro Dalcin   jac->tol                = PETSC_DEFAULT;
76116d9e3a6SLisandro Dalcin   jac->printstatistics    = PetscLogPrintInfo;
76216d9e3a6SLisandro Dalcin 
76316d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr);
76416d9e3a6SLisandro Dalcin   if (flag) {
76530e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParCSRPilutCreate,(jac->comm_hypre,&jac->hsolver));
76616d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
76716d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Pilut;
76816d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParCSRPilutDestroy;
76916d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParCSRPilutSetup;
77016d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParCSRPilutSolve;
77116d9e3a6SLisandro Dalcin     jac->factorrowsize      = PETSC_DEFAULT;
77216d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
77316d9e3a6SLisandro Dalcin   }
77416d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr);
77516d9e3a6SLisandro Dalcin   if (flag) {
77630e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsCreate,(jac->comm_hypre,&jac->hsolver));
77716d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
77816d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_ParaSails;
77916d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParaSailsDestroy;
78016d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParaSailsSetup;
78116d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParaSailsSolve;
78216d9e3a6SLisandro Dalcin     /* initialize */
78316d9e3a6SLisandro Dalcin     jac->nlevels     = 1;
78416d9e3a6SLisandro Dalcin     jac->threshhold  = .1;
78516d9e3a6SLisandro Dalcin     jac->filter      = .1;
78616d9e3a6SLisandro Dalcin     jac->loadbal     = 0;
78716d9e3a6SLisandro Dalcin     if (PetscLogPrintInfo) {
78816d9e3a6SLisandro Dalcin       jac->logging     = (int) PETSC_TRUE;
78916d9e3a6SLisandro Dalcin     } else {
79016d9e3a6SLisandro Dalcin       jac->logging     = (int) PETSC_FALSE;
79116d9e3a6SLisandro Dalcin     }
79216d9e3a6SLisandro Dalcin     jac->ruse = (int) PETSC_FALSE;
79316d9e3a6SLisandro Dalcin     jac->symt   = 0;
79430e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels));
79530e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter));
79630e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal));
79730e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging));
79830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse));
79930e6f737SJed Brown     PetscStackCallHypre(0,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;
85130e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype));
85230e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels));
85330e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter));
85430e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol));
85530e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor));
85630e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold));
85730e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum));
85830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype));
85930e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype));
86030e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder));
86130e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype));
86230e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl));
86330e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax));
86430e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths));
86530e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, jac->relaxtype[0]));  /*defaults coarse to 9*/
86630e6f737SJed Brown     PetscStackCallHypre(0,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);
87016d9e3a6SLisandro Dalcin   jac->hypre_type = PETSC_NULL;
87165e19b50SBarry Smith   SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
87216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
87316d9e3a6SLisandro Dalcin }
87416d9e3a6SLisandro Dalcin EXTERN_C_END
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 EXTERN_C_BEGIN
100116d9e3a6SLisandro Dalcin #undef __FUNCT__
100216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCCreate_HYPRE"
10037087cfbeSBarry Smith PetscErrorCode  PCCreate_HYPRE(PC pc)
100416d9e3a6SLisandro Dalcin {
100516d9e3a6SLisandro Dalcin   PC_HYPRE       *jac;
100616d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
100716d9e3a6SLisandro Dalcin 
100816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
100938f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_HYPRE,&jac);CHKERRQ(ierr);
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;
101616d9e3a6SLisandro Dalcin   jac->hypre_type          = PETSC_NULL;
101716d9e3a6SLisandro Dalcin   /* duplicate communicator for hypre */
10187adad957SLisandro Dalcin   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(jac->comm_hypre));CHKERRQ(ierr);
101916d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr);
102016d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr);
102116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
102216d9e3a6SLisandro Dalcin }
102316d9e3a6SLisandro Dalcin EXTERN_C_END
1024ebc551c0SBarry Smith 
1025f91d8e95SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------*/
1026f91d8e95SBarry Smith 
1027b862ddfaSBarry Smith /* this include is needed ONLY to allow access to the private data inside the Mat object specific to hypre */
1028*b45d2f2cSJed Brown #include <petsc-private/matimpl.h>
1029ebc551c0SBarry Smith 
1030ebc551c0SBarry Smith typedef struct {
103168326731SBarry Smith   MPI_Comm            hcomm;       /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */
1032f91d8e95SBarry Smith   HYPRE_StructSolver  hsolver;
10339e5bc791SBarry Smith 
10349e5bc791SBarry Smith   /* keep copy of PFMG options used so may view them */
10354ddd07fcSJed Brown   PetscInt            its;
10369e5bc791SBarry Smith   double              tol;
10374ddd07fcSJed Brown   PetscInt            relax_type;
10384ddd07fcSJed Brown   PetscInt            rap_type;
10394ddd07fcSJed Brown   PetscInt            num_pre_relax,num_post_relax;
10404ddd07fcSJed Brown   PetscInt            max_levels;
1041ebc551c0SBarry Smith } PC_PFMG;
1042ebc551c0SBarry Smith 
1043ebc551c0SBarry Smith #undef __FUNCT__
1044ebc551c0SBarry Smith #define __FUNCT__ "PCDestroy_PFMG"
1045ebc551c0SBarry Smith PetscErrorCode PCDestroy_PFMG(PC pc)
1046ebc551c0SBarry Smith {
1047ebc551c0SBarry Smith   PetscErrorCode ierr;
1048f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1049ebc551c0SBarry Smith 
1050ebc551c0SBarry Smith   PetscFunctionBegin;
105130e6f737SJed Brown   if (ex->hsolver) {PetscStackCallHypre(0,HYPRE_StructPFMGDestroy,(ex->hsolver));}
1052f91d8e95SBarry Smith   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1053c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1054ebc551c0SBarry Smith   PetscFunctionReturn(0);
1055ebc551c0SBarry Smith }
1056ebc551c0SBarry Smith 
10579e5bc791SBarry Smith static const char *PFMGRelaxType[]   = {"Jacobi","Weighted-Jacobi","symmetric-Red/Black-Gauss-Seidel","Red/Black-Gauss-Seidel"};
10589e5bc791SBarry Smith static const char *PFMGRAPType[]   = {"Galerkin","non-Galerkin"};
10599e5bc791SBarry Smith 
1060ebc551c0SBarry Smith #undef __FUNCT__
1061ebc551c0SBarry Smith #define __FUNCT__ "PCView_PFMG"
1062ebc551c0SBarry Smith PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer)
1063ebc551c0SBarry Smith {
1064ebc551c0SBarry Smith   PetscErrorCode ierr;
1065ace3abfcSBarry Smith   PetscBool      iascii;
1066f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1067ebc551c0SBarry Smith 
1068ebc551c0SBarry Smith   PetscFunctionBegin;
10692692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10709e5bc791SBarry Smith   if (iascii) {
10719e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG preconditioning\n");CHKERRQ(ierr);
10729e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
10739e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
10749e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
10759e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: RAP type %s\n",PFMGRAPType[ex->rap_type]);CHKERRQ(ierr);
10769e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
10773b46a515SGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max levels %d\n",ex->max_levels);CHKERRQ(ierr);
10789e5bc791SBarry Smith   }
1079ebc551c0SBarry Smith   PetscFunctionReturn(0);
1080ebc551c0SBarry Smith }
1081ebc551c0SBarry Smith 
10829e5bc791SBarry Smith 
1083ebc551c0SBarry Smith #undef __FUNCT__
1084ebc551c0SBarry Smith #define __FUNCT__ "PCSetFromOptions_PFMG"
1085ebc551c0SBarry Smith PetscErrorCode PCSetFromOptions_PFMG(PC pc)
1086ebc551c0SBarry Smith {
1087ebc551c0SBarry Smith   PetscErrorCode ierr;
1088f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1089ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
1090ebc551c0SBarry Smith 
1091ebc551c0SBarry Smith   PetscFunctionBegin;
1092ebc551c0SBarry Smith   ierr = PetscOptionsHead("PFMG options");CHKERRQ(ierr);
1093acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
109468326731SBarry Smith   if (flg) {
1095a0324ebeSBarry Smith     int level=3;
109630e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_StructPFMGSetPrintLevel,(ex->hsolver,level));
109768326731SBarry Smith   }
10989e5bc791SBarry Smith   ierr = PetscOptionsInt("-pc_pfmg_its","Number of iterations of PFMG to use as preconditioner","HYPRE_StructPFMGSetMaxIter",ex->its,&ex->its,PETSC_NULL);CHKERRQ(ierr);
109930e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetMaxIter,(ex->hsolver,ex->its));
11009e5bc791SBarry Smith   ierr = PetscOptionsInt("-pc_pfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_StructPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,PETSC_NULL);CHKERRQ(ierr);
110130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetNumPreRelax,(ex->hsolver,ex->num_pre_relax));
11029e5bc791SBarry Smith   ierr = PetscOptionsInt("-pc_pfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_StructPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,PETSC_NULL);CHKERRQ(ierr);
110330e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetNumPostRelax,(ex->hsolver,ex->num_post_relax));
11049e5bc791SBarry Smith 
11053b46a515SGlenn Hammond   ierr = PetscOptionsInt("-pc_pfmg_max_levels","Max Levels for MG hierarchy","HYPRE_StructPFMGSetMaxLevels",ex->max_levels,&ex->max_levels,PETSC_NULL);CHKERRQ(ierr);
110630e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetMaxLevels,(ex->hsolver,ex->max_levels));
11073b46a515SGlenn Hammond 
11089e5bc791SBarry Smith   ierr = PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr);
110930e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetTol,(ex->hsolver,ex->tol));
11109e5bc791SBarry Smith   ierr = PetscOptionsEList("-pc_pfmg_relax_type","Relax type for the up and down cycles","HYPRE_StructPFMGSetRelaxType",PFMGRelaxType,4,PFMGRelaxType[ex->relax_type],&ex->relax_type,PETSC_NULL);CHKERRQ(ierr);
111130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetRelaxType,(ex->hsolver, ex->relax_type));
11129e5bc791SBarry Smith   ierr = PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,2,PFMGRAPType[ex->rap_type],&ex->rap_type,PETSC_NULL);CHKERRQ(ierr);
111330e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetRAPType,(ex->hsolver, ex->rap_type));
1114ebc551c0SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1115ebc551c0SBarry Smith   PetscFunctionReturn(0);
1116ebc551c0SBarry Smith }
1117ebc551c0SBarry Smith 
1118f91d8e95SBarry Smith #undef __FUNCT__
1119f91d8e95SBarry Smith #define __FUNCT__ "PCApply_PFMG"
1120f91d8e95SBarry Smith PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y)
1121f91d8e95SBarry Smith {
1122f91d8e95SBarry Smith   PetscErrorCode  ierr;
1123f91d8e95SBarry Smith   PC_PFMG         *ex = (PC_PFMG*) pc->data;
1124f91d8e95SBarry Smith   PetscScalar     *xx,*yy;
11254ddd07fcSJed Brown   PetscInt        ilower[3],iupper[3];
112668326731SBarry Smith   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
1127f91d8e95SBarry Smith 
1128f91d8e95SBarry Smith   PetscFunctionBegin;
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 */
113530e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructVectorSetConstantValues,(mx->hb,0.0));
1136f91d8e95SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
113730e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructVectorSetBoxValues,(mx->hb,ilower,iupper,xx));
1138f91d8e95SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
113930e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructVectorAssemble,(mx->hb));
114030e6f737SJed Brown   PetscStackCallHypre(0,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);
114430e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructVectorGetBoxValues,(mx->hx,ilower,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;
115830e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetMaxIter,(jac->hsolver,its*jac->its));
115930e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetTol,(jac->hsolver,rtol));
11609e5bc791SBarry Smith 
11619e5bc791SBarry Smith   ierr = PCApply_PFMG(pc,b,y);CHKERRQ(ierr);
116230e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGGetNumIterations,(jac->hsolver,&oits));
11639e5bc791SBarry Smith   *outits = oits;
11649e5bc791SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
11659e5bc791SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
116630e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetTol,(jac->hsolver,jac->tol));
116730e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetMaxIter,(jac->hsolver,jac->its));
11689e5bc791SBarry Smith   PetscFunctionReturn(0);
11699e5bc791SBarry Smith }
11709e5bc791SBarry Smith 
11719e5bc791SBarry Smith 
11723a32d3dbSGlenn Hammond #undef __FUNCT__
11733a32d3dbSGlenn Hammond #define __FUNCT__ "PCSetUp_PFMG"
11743a32d3dbSGlenn Hammond PetscErrorCode PCSetUp_PFMG(PC pc)
11753a32d3dbSGlenn Hammond {
11763a32d3dbSGlenn Hammond   PetscErrorCode  ierr;
11773a32d3dbSGlenn Hammond   PC_PFMG         *ex = (PC_PFMG*) pc->data;
11783a32d3dbSGlenn Hammond   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
1179ace3abfcSBarry Smith   PetscBool       flg;
11803a32d3dbSGlenn Hammond 
11813a32d3dbSGlenn Hammond   PetscFunctionBegin;
11823a32d3dbSGlenn Hammond   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATHYPRESTRUCT,&flg);CHKERRQ(ierr);
1183e7e72b3dSBarry Smith   if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESTRUCT with this preconditioner");
11843a32d3dbSGlenn Hammond 
11853a32d3dbSGlenn Hammond   /* create the hypre solver object and set its information */
11863a32d3dbSGlenn Hammond   if (ex->hsolver) {
118730e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_StructPFMGDestroy,(ex->hsolver));
11883a32d3dbSGlenn Hammond   }
118930e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
11903a32d3dbSGlenn Hammond   ierr = PCSetFromOptions_PFMG(pc);CHKERRQ(ierr);
119130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetup,(ex->hsolver,mx->hmat,mx->hb,mx->hx));
119230e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGSetZeroGuess,(ex->hsolver));
11933a32d3dbSGlenn Hammond   PetscFunctionReturn(0);
11943a32d3dbSGlenn Hammond }
11953a32d3dbSGlenn Hammond 
1196ebc551c0SBarry Smith 
1197ebc551c0SBarry Smith /*MC
1198ebc551c0SBarry Smith      PCPFMG - the hypre PFMG multigrid solver
1199ebc551c0SBarry Smith 
1200ebc551c0SBarry Smith    Level: advanced
1201ebc551c0SBarry Smith 
12029e5bc791SBarry Smith    Options Database:
12039e5bc791SBarry Smith + -pc_pfmg_its <its> number of iterations of PFMG to use as preconditioner
12049e5bc791SBarry Smith . -pc_pfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
12059e5bc791SBarry Smith . -pc_pfmg_num_post_relax <steps> number of smoothing steps after coarse grid
12069e5bc791SBarry Smith . -pc_pfmg_tol <tol> tolerance of PFMG
12079e5bc791SBarry 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
12089e5bc791SBarry Smith - -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin
1209f91d8e95SBarry Smith 
12109e5bc791SBarry Smith    Notes:  This is for CELL-centered descretizations
12119e5bc791SBarry Smith 
12128e395302SJed Brown            This must be used with the MATHYPRESTRUCT matrix type.
1213aa219208SBarry Smith            This is less general than in hypre, it supports only one block per process defined by a PETSc DMDA.
12149e5bc791SBarry Smith 
12159e5bc791SBarry Smith .seealso:  PCMG, MATHYPRESTRUCT
1216ebc551c0SBarry Smith M*/
1217ebc551c0SBarry Smith 
1218ebc551c0SBarry Smith EXTERN_C_BEGIN
1219ebc551c0SBarry Smith #undef __FUNCT__
1220ebc551c0SBarry Smith #define __FUNCT__ "PCCreate_PFMG"
12217087cfbeSBarry Smith PetscErrorCode  PCCreate_PFMG(PC pc)
1222ebc551c0SBarry Smith {
1223ebc551c0SBarry Smith   PetscErrorCode ierr;
1224ebc551c0SBarry Smith   PC_PFMG        *ex;
1225ebc551c0SBarry Smith 
1226ebc551c0SBarry Smith   PetscFunctionBegin;
1227ebc551c0SBarry Smith   ierr = PetscNew(PC_PFMG,&ex);CHKERRQ(ierr);\
122868326731SBarry Smith   pc->data = ex;
1229ebc551c0SBarry Smith 
12309e5bc791SBarry Smith   ex->its            = 1;
12319e5bc791SBarry Smith   ex->tol            = 1.e-8;
12329e5bc791SBarry Smith   ex->relax_type     = 1;
12339e5bc791SBarry Smith   ex->rap_type       = 0;
12349e5bc791SBarry Smith   ex->num_pre_relax  = 1;
12359e5bc791SBarry Smith   ex->num_post_relax = 1;
12363b46a515SGlenn Hammond   ex->max_levels     = 0;
12379e5bc791SBarry Smith 
1238ebc551c0SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_PFMG;
1239ebc551c0SBarry Smith   pc->ops->view            = PCView_PFMG;
1240ebc551c0SBarry Smith   pc->ops->destroy         = PCDestroy_PFMG;
1241f91d8e95SBarry Smith   pc->ops->apply           = PCApply_PFMG;
12429e5bc791SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_PFMG;
124368326731SBarry Smith   pc->ops->setup           = PCSetUp_PFMG;
1244f91d8e95SBarry Smith   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr);
124530e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver));
1246ebc551c0SBarry Smith   PetscFunctionReturn(0);
1247ebc551c0SBarry Smith }
1248ebc551c0SBarry Smith EXTERN_C_END
1249d851a50bSGlenn Hammond 
1250325fc9f4SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------------------------*/
1251325fc9f4SBarry Smith 
1252d851a50bSGlenn Hammond /* we know we are working with a HYPRE_SStructMatrix */
1253d851a50bSGlenn Hammond typedef struct {
1254d851a50bSGlenn Hammond   MPI_Comm            hcomm;       /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */
1255d851a50bSGlenn Hammond   HYPRE_SStructSolver  ss_solver;
1256d851a50bSGlenn Hammond 
1257d851a50bSGlenn Hammond   /* keep copy of SYSPFMG options used so may view them */
12584ddd07fcSJed Brown   PetscInt            its;
1259d851a50bSGlenn Hammond   double              tol;
12604ddd07fcSJed Brown   PetscInt            relax_type;
12614ddd07fcSJed Brown   PetscInt            num_pre_relax,num_post_relax;
1262d851a50bSGlenn Hammond } PC_SysPFMG;
1263d851a50bSGlenn Hammond 
1264d851a50bSGlenn Hammond #undef __FUNCT__
1265d851a50bSGlenn Hammond #define __FUNCT__ "PCDestroy_SysPFMG"
1266d851a50bSGlenn Hammond PetscErrorCode PCDestroy_SysPFMG(PC pc)
1267d851a50bSGlenn Hammond {
1268d851a50bSGlenn Hammond   PetscErrorCode ierr;
1269d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1270d851a50bSGlenn Hammond 
1271d851a50bSGlenn Hammond   PetscFunctionBegin;
127230e6f737SJed Brown   if (ex->ss_solver) {PetscStackCallHypre(0,HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));}
1273d851a50bSGlenn Hammond   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1274c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
1275d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1276d851a50bSGlenn Hammond }
1277d851a50bSGlenn Hammond 
1278d851a50bSGlenn Hammond static const char *SysPFMGRelaxType[]   = {"Weighted-Jacobi","Red/Black-Gauss-Seidel"};
1279d851a50bSGlenn Hammond 
1280d851a50bSGlenn Hammond #undef __FUNCT__
1281d851a50bSGlenn Hammond #define __FUNCT__ "PCView_SysPFMG"
1282d851a50bSGlenn Hammond PetscErrorCode PCView_SysPFMG(PC pc,PetscViewer viewer)
1283d851a50bSGlenn Hammond {
1284d851a50bSGlenn Hammond   PetscErrorCode ierr;
1285ace3abfcSBarry Smith   PetscBool      iascii;
1286d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1287d851a50bSGlenn Hammond 
1288d851a50bSGlenn Hammond   PetscFunctionBegin;
12892692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1290d851a50bSGlenn Hammond   if (iascii) {
1291d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG preconditioning\n");CHKERRQ(ierr);
1292d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
1293d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
1294d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
1295d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
1296d851a50bSGlenn Hammond   }
1297d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1298d851a50bSGlenn Hammond }
1299d851a50bSGlenn Hammond 
1300d851a50bSGlenn Hammond 
1301d851a50bSGlenn Hammond #undef __FUNCT__
1302d851a50bSGlenn Hammond #define __FUNCT__ "PCSetFromOptions_SysPFMG"
1303d851a50bSGlenn Hammond PetscErrorCode PCSetFromOptions_SysPFMG(PC pc)
1304d851a50bSGlenn Hammond {
1305d851a50bSGlenn Hammond   PetscErrorCode ierr;
1306d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1307ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
1308d851a50bSGlenn Hammond 
1309d851a50bSGlenn Hammond   PetscFunctionBegin;
1310d851a50bSGlenn Hammond   ierr = PetscOptionsHead("SysPFMG options");CHKERRQ(ierr);
1311acfcf0e5SJed Brown   ierr = PetscOptionsBool("-pc_syspfmg_print_statistics","Print statistics","HYPRE_SStructSysPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1312d851a50bSGlenn Hammond   if (flg) {
1313d851a50bSGlenn Hammond     int level=3;
131430e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetPrintLevel,(ex->ss_solver,level));
1315d851a50bSGlenn Hammond   }
1316d851a50bSGlenn Hammond   ierr = PetscOptionsInt("-pc_syspfmg_its","Number of iterations of SysPFMG to use as preconditioner","HYPRE_SStructSysPFMGSetMaxIter",ex->its,&ex->its,PETSC_NULL);CHKERRQ(ierr);
131730e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetMaxIter,(ex->ss_solver,ex->its));
1318d851a50bSGlenn Hammond   ierr = PetscOptionsInt("-pc_syspfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_SStructSysPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,PETSC_NULL);CHKERRQ(ierr);
131930e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetNumPreRelax,(ex->ss_solver,ex->num_pre_relax));
1320d851a50bSGlenn Hammond   ierr = PetscOptionsInt("-pc_syspfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_SStructSysPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,PETSC_NULL);CHKERRQ(ierr);
132130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetNumPostRelax,(ex->ss_solver,ex->num_post_relax));
1322d851a50bSGlenn Hammond 
1323d851a50bSGlenn Hammond   ierr = PetscOptionsReal("-pc_syspfmg_tol","Tolerance of SysPFMG","HYPRE_SStructSysPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr);
132430e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetTol,(ex->ss_solver,ex->tol));
1325d851a50bSGlenn Hammond   ierr = PetscOptionsEList("-pc_syspfmg_relax_type","Relax type for the up and down cycles","HYPRE_SStructSysPFMGSetRelaxType",SysPFMGRelaxType,4,SysPFMGRelaxType[ex->relax_type],&ex->relax_type,PETSC_NULL);CHKERRQ(ierr);
132630e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetRelaxType,(ex->ss_solver, ex->relax_type));
1327d851a50bSGlenn Hammond   ierr = PetscOptionsTail();CHKERRQ(ierr);
1328d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1329d851a50bSGlenn Hammond }
1330d851a50bSGlenn Hammond 
1331d851a50bSGlenn Hammond #undef __FUNCT__
1332d851a50bSGlenn Hammond #define __FUNCT__ "PCApply_SysPFMG"
1333d851a50bSGlenn Hammond PetscErrorCode PCApply_SysPFMG(PC pc,Vec x,Vec y)
1334d851a50bSGlenn Hammond {
1335d851a50bSGlenn Hammond   PetscErrorCode    ierr;
1336d851a50bSGlenn Hammond   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
1337d851a50bSGlenn Hammond   PetscScalar      *xx,*yy;
13384ddd07fcSJed Brown   PetscInt          ilower[3],iupper[3];
1339d851a50bSGlenn Hammond   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct *)(pc->pmat->data);
13404ddd07fcSJed Brown   PetscInt          ordering= mx->dofs_order;
13414ddd07fcSJed Brown   PetscInt          nvars= mx->nvars;
13424ddd07fcSJed Brown   PetscInt          part= 0;
13434ddd07fcSJed Brown   PetscInt          size;
13444ddd07fcSJed Brown   PetscInt          i;
1345d851a50bSGlenn Hammond 
1346d851a50bSGlenn Hammond   PetscFunctionBegin;
1347aa219208SBarry Smith   ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1348d851a50bSGlenn Hammond   iupper[0] += ilower[0] - 1;
1349d851a50bSGlenn Hammond   iupper[1] += ilower[1] - 1;
1350d851a50bSGlenn Hammond   iupper[2] += ilower[2] - 1;
1351d851a50bSGlenn Hammond 
1352d851a50bSGlenn Hammond   size= 1;
1353d851a50bSGlenn Hammond   for (i= 0; i< 3; i++) {
1354d851a50bSGlenn Hammond      size*= (iupper[i]-ilower[i]+1);
1355d851a50bSGlenn Hammond   }
1356d851a50bSGlenn Hammond   /* copy x values over to hypre for variable ordering */
1357d851a50bSGlenn Hammond   if (ordering) {
135830e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
1359d851a50bSGlenn Hammond      ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1360d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
136130e6f737SJed Brown        PetscStackCallHypre(0,HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,ilower,iupper,i,xx+(size*i)));
1362d851a50bSGlenn Hammond      }
1363d851a50bSGlenn Hammond      ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
136430e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_SStructVectorAssemble,(mx->ss_b));
136530e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
136630e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1367d851a50bSGlenn Hammond 
1368d851a50bSGlenn Hammond      /* copy solution values back to PETSc */
1369d851a50bSGlenn Hammond      ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1370d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
137130e6f737SJed Brown        PetscStackCallHypre(0,HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,ilower,iupper,i,yy+(size*i)));
1372d851a50bSGlenn Hammond      }
1373d851a50bSGlenn Hammond      ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1374a65764d7SBarry Smith   } else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
1375d851a50bSGlenn Hammond      PetscScalar     *z;
13764ddd07fcSJed Brown      PetscInt         j, k;
1377d851a50bSGlenn Hammond 
1378d851a50bSGlenn Hammond      ierr = PetscMalloc(nvars*size*sizeof(PetscScalar),&z);CHKERRQ(ierr);
137930e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
1380d851a50bSGlenn Hammond      ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1381d851a50bSGlenn Hammond 
1382d851a50bSGlenn Hammond      /* transform nodal to hypre's variable ordering for sys_pfmg */
1383d851a50bSGlenn Hammond      for (i= 0; i< size; i++) {
1384d851a50bSGlenn Hammond         k= i*nvars;
1385d851a50bSGlenn Hammond         for (j= 0; j< nvars; j++) {
1386d851a50bSGlenn Hammond            z[j*size+i]= xx[k+j];
1387d851a50bSGlenn Hammond         }
1388d851a50bSGlenn Hammond      }
1389d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
139030e6f737SJed Brown        PetscStackCallHypre(0,HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,ilower,iupper,i,z+(size*i)));
1391d851a50bSGlenn Hammond      }
1392d851a50bSGlenn Hammond      ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
139330e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_SStructVectorAssemble,(mx->ss_b));
139430e6f737SJed Brown      PetscStackCallHypre(0,HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1395d851a50bSGlenn Hammond 
1396d851a50bSGlenn Hammond      /* copy solution values back to PETSc */
1397d851a50bSGlenn Hammond      ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1398d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
139930e6f737SJed Brown        PetscStackCallHypre(0,HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,ilower,iupper,i,z+(size*i)));
1400d851a50bSGlenn Hammond      }
1401d851a50bSGlenn Hammond      /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
1402d851a50bSGlenn Hammond      for (i= 0; i< size; i++) {
1403d851a50bSGlenn Hammond         k= i*nvars;
1404d851a50bSGlenn Hammond         for (j= 0; j< nvars; j++) {
1405d851a50bSGlenn Hammond            yy[k+j]= z[j*size+i];
1406d851a50bSGlenn Hammond         }
1407d851a50bSGlenn Hammond      }
1408d851a50bSGlenn Hammond      ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1409d851a50bSGlenn Hammond      ierr = PetscFree(z);CHKERRQ(ierr);
1410d851a50bSGlenn Hammond   }
1411d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1412d851a50bSGlenn Hammond }
1413d851a50bSGlenn Hammond 
1414d851a50bSGlenn Hammond #undef __FUNCT__
1415d851a50bSGlenn Hammond #define __FUNCT__ "PCApplyRichardson_SysPFMG"
1416ace3abfcSBarry 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)
1417d851a50bSGlenn Hammond {
1418d851a50bSGlenn Hammond   PC_SysPFMG    *jac = (PC_SysPFMG*)pc->data;
1419d851a50bSGlenn Hammond   PetscErrorCode ierr;
14204ddd07fcSJed Brown   PetscInt       oits;
1421d851a50bSGlenn Hammond 
1422d851a50bSGlenn Hammond   PetscFunctionBegin;
142330e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,its*jac->its));
142430e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,rtol));
1425d851a50bSGlenn Hammond   ierr = PCApply_SysPFMG(pc,b,y);CHKERRQ(ierr);
142630e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGGetNumIterations,(jac->ss_solver,&oits));
1427d851a50bSGlenn Hammond   *outits = oits;
1428d851a50bSGlenn Hammond   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
1429d851a50bSGlenn Hammond   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
143030e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,jac->tol));
143130e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,jac->its));
1432d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1433d851a50bSGlenn Hammond }
1434d851a50bSGlenn Hammond 
1435d851a50bSGlenn Hammond 
1436d851a50bSGlenn Hammond #undef __FUNCT__
1437d851a50bSGlenn Hammond #define __FUNCT__ "PCSetUp_SysPFMG"
1438d851a50bSGlenn Hammond PetscErrorCode PCSetUp_SysPFMG(PC pc)
1439d851a50bSGlenn Hammond {
1440d851a50bSGlenn Hammond   PetscErrorCode    ierr;
1441d851a50bSGlenn Hammond   PC_SysPFMG        *ex = (PC_SysPFMG*) pc->data;
1442d851a50bSGlenn Hammond   Mat_HYPRESStruct  *mx = (Mat_HYPRESStruct *)(pc->pmat->data);
1443ace3abfcSBarry Smith   PetscBool         flg;
1444d851a50bSGlenn Hammond 
1445d851a50bSGlenn Hammond   PetscFunctionBegin;
1446d851a50bSGlenn Hammond   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATHYPRESSTRUCT,&flg);CHKERRQ(ierr);
1447e7e72b3dSBarry Smith   if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESSTRUCT with this preconditioner");
1448d851a50bSGlenn Hammond 
1449d851a50bSGlenn Hammond   /* create the hypre sstruct solver object and set its information */
1450d851a50bSGlenn Hammond   if (ex->ss_solver) {
145130e6f737SJed Brown     PetscStackCallHypre(0,HYPRE_SStructSysPFMGDestroy,(ex->ss_solver));
1452d851a50bSGlenn Hammond   }
145330e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
1454d851a50bSGlenn Hammond   ierr = PCSetFromOptions_SysPFMG(pc);CHKERRQ(ierr);
145530e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetZeroGuess,(ex->ss_solver));
145630e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGSetup,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1457d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1458d851a50bSGlenn Hammond }
1459d851a50bSGlenn Hammond 
1460d851a50bSGlenn Hammond 
1461d851a50bSGlenn Hammond /*MC
1462d851a50bSGlenn Hammond      PCSysPFMG - the hypre SysPFMG multigrid solver
1463d851a50bSGlenn Hammond 
1464d851a50bSGlenn Hammond    Level: advanced
1465d851a50bSGlenn Hammond 
1466d851a50bSGlenn Hammond    Options Database:
1467d851a50bSGlenn Hammond + -pc_syspfmg_its <its> number of iterations of SysPFMG to use as preconditioner
1468d851a50bSGlenn Hammond . -pc_syspfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
1469d851a50bSGlenn Hammond . -pc_syspfmg_num_post_relax <steps> number of smoothing steps after coarse grid
1470d851a50bSGlenn Hammond . -pc_syspfmg_tol <tol> tolerance of SysPFMG
1471d851a50bSGlenn Hammond . -pc_syspfmg_relax_type -relaxation type for the up and down cycles, one of Weighted-Jacobi,Red/Black-Gauss-Seidel
1472d851a50bSGlenn Hammond 
1473d851a50bSGlenn Hammond    Notes:  This is for CELL-centered descretizations
1474d851a50bSGlenn Hammond 
1475f6680f47SSatish Balay            This must be used with the MATHYPRESSTRUCT matrix type.
1476aa219208SBarry Smith            This is less general than in hypre, it supports only one part, and one block per process defined by a PETSc DMDA.
1477d851a50bSGlenn Hammond            Also, only cell-centered variables.
1478d851a50bSGlenn Hammond 
1479d851a50bSGlenn Hammond .seealso:  PCMG, MATHYPRESSTRUCT
1480d851a50bSGlenn Hammond M*/
1481d851a50bSGlenn Hammond 
1482d851a50bSGlenn Hammond EXTERN_C_BEGIN
1483d851a50bSGlenn Hammond #undef __FUNCT__
1484d851a50bSGlenn Hammond #define __FUNCT__ "PCCreate_SysPFMG"
14857087cfbeSBarry Smith PetscErrorCode  PCCreate_SysPFMG(PC pc)
1486d851a50bSGlenn Hammond {
1487d851a50bSGlenn Hammond   PetscErrorCode     ierr;
1488d851a50bSGlenn Hammond   PC_SysPFMG        *ex;
1489d851a50bSGlenn Hammond 
1490d851a50bSGlenn Hammond   PetscFunctionBegin;
1491d851a50bSGlenn Hammond   ierr = PetscNew(PC_SysPFMG,&ex);CHKERRQ(ierr);\
1492d851a50bSGlenn Hammond   pc->data = ex;
1493d851a50bSGlenn Hammond 
1494d851a50bSGlenn Hammond   ex->its            = 1;
1495d851a50bSGlenn Hammond   ex->tol            = 1.e-8;
1496d851a50bSGlenn Hammond   ex->relax_type     = 1;
1497d851a50bSGlenn Hammond   ex->num_pre_relax  = 1;
1498d851a50bSGlenn Hammond   ex->num_post_relax = 1;
1499d851a50bSGlenn Hammond 
1500d851a50bSGlenn Hammond   pc->ops->setfromoptions  = PCSetFromOptions_SysPFMG;
1501d851a50bSGlenn Hammond   pc->ops->view            = PCView_SysPFMG;
1502d851a50bSGlenn Hammond   pc->ops->destroy         = PCDestroy_SysPFMG;
1503d851a50bSGlenn Hammond   pc->ops->apply           = PCApply_SysPFMG;
1504d851a50bSGlenn Hammond   pc->ops->applyrichardson = PCApplyRichardson_SysPFMG;
1505d851a50bSGlenn Hammond   pc->ops->setup           = PCSetUp_SysPFMG;
1506d851a50bSGlenn Hammond   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr);
150730e6f737SJed Brown   PetscStackCallHypre(0,HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver));
1508d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1509d851a50bSGlenn Hammond }
1510d851a50bSGlenn Hammond EXTERN_C_END
1511