xref: /petsc/src/ksp/pc/impls/hypre/hypre.c (revision a65764d7969f078334d4978dbb46d5ea9f0bcbcf)
116d9e3a6SLisandro Dalcin #define PETSCKSP_DLL
216d9e3a6SLisandro Dalcin 
316d9e3a6SLisandro Dalcin /*
416d9e3a6SLisandro Dalcin    Provides an interface to the LLNL package hypre
516d9e3a6SLisandro Dalcin */
60f1074feSSatish Balay 
70f1074feSSatish Balay /* Must use hypre 2.0.0 or more recent. */
80f1074feSSatish Balay 
916d9e3a6SLisandro Dalcin #include "private/pcimpl.h"          /*I "petscpc.h" I*/
109f73f8ecSBarry Smith #include "../src/dm/da/utils/mhyp.h"
1116d9e3a6SLisandro Dalcin 
1216d9e3a6SLisandro Dalcin /*
1316d9e3a6SLisandro Dalcin    Private context (data structure) for the  preconditioner.
1416d9e3a6SLisandro Dalcin */
1516d9e3a6SLisandro Dalcin typedef struct {
1616d9e3a6SLisandro Dalcin   HYPRE_Solver       hsolver;
1716d9e3a6SLisandro Dalcin   HYPRE_IJMatrix     ij;
1816d9e3a6SLisandro Dalcin   HYPRE_IJVector     b,x;
1916d9e3a6SLisandro Dalcin 
2016d9e3a6SLisandro Dalcin   PetscErrorCode     (*destroy)(HYPRE_Solver);
2116d9e3a6SLisandro Dalcin   PetscErrorCode     (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
2216d9e3a6SLisandro Dalcin   PetscErrorCode     (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector);
2316d9e3a6SLisandro Dalcin 
2416d9e3a6SLisandro Dalcin   MPI_Comm           comm_hypre;
2516d9e3a6SLisandro Dalcin   char              *hypre_type;
2616d9e3a6SLisandro Dalcin 
2716d9e3a6SLisandro Dalcin   /* options for Pilut and BoomerAMG*/
2816d9e3a6SLisandro Dalcin   int                maxiter;
2916d9e3a6SLisandro Dalcin   double             tol;
3016d9e3a6SLisandro Dalcin 
3116d9e3a6SLisandro Dalcin   /* options for Pilut */
3216d9e3a6SLisandro Dalcin   int                factorrowsize;
3316d9e3a6SLisandro Dalcin 
3416d9e3a6SLisandro Dalcin   /* options for ParaSails */
3516d9e3a6SLisandro Dalcin   int                nlevels;
3616d9e3a6SLisandro Dalcin   double             threshhold;
3716d9e3a6SLisandro Dalcin   double             filter;
3816d9e3a6SLisandro Dalcin   int                sym;
3916d9e3a6SLisandro Dalcin   double             loadbal;
4016d9e3a6SLisandro Dalcin   int                logging;
4116d9e3a6SLisandro Dalcin   int                ruse;
4216d9e3a6SLisandro Dalcin   int                symt;
4316d9e3a6SLisandro Dalcin 
4416d9e3a6SLisandro Dalcin   /* options for Euclid */
4516d9e3a6SLisandro Dalcin   PetscTruth         bjilu;
4616d9e3a6SLisandro Dalcin   int                levels;
4716d9e3a6SLisandro Dalcin 
4816d9e3a6SLisandro Dalcin   /* options for Euclid and BoomerAMG */
4916d9e3a6SLisandro Dalcin   PetscTruth         printstatistics;
5016d9e3a6SLisandro Dalcin 
5116d9e3a6SLisandro Dalcin   /* options for BoomerAMG */
5216d9e3a6SLisandro Dalcin   int                cycletype;
5316d9e3a6SLisandro Dalcin   int                maxlevels;
5416d9e3a6SLisandro Dalcin   double             strongthreshold;
5516d9e3a6SLisandro Dalcin   double             maxrowsum;
560f1074feSSatish Balay   int                gridsweeps[3];
5716d9e3a6SLisandro Dalcin   int                coarsentype;
5816d9e3a6SLisandro Dalcin   int                measuretype;
590f1074feSSatish Balay   int                relaxtype[3];
6016d9e3a6SLisandro Dalcin   double             relaxweight;
6116d9e3a6SLisandro Dalcin   double             outerrelaxweight;
6216d9e3a6SLisandro Dalcin   int                relaxorder;
6316d9e3a6SLisandro Dalcin   double             truncfactor;
6416d9e3a6SLisandro Dalcin   PetscTruth         applyrichardson;
650f1074feSSatish Balay   int                pmax;
660f1074feSSatish Balay   int                interptype;
670f1074feSSatish Balay   int                agg_nl;
680f1074feSSatish Balay   int                agg_num_paths;
698f87f92bSBarry Smith   int                nodal_coarsen;
708f87f92bSBarry Smith   PetscTruth         nodal_relax;
718f87f92bSBarry Smith   int                nodal_relax_levels;
7216d9e3a6SLisandro Dalcin } PC_HYPRE;
7316d9e3a6SLisandro Dalcin 
7416d9e3a6SLisandro Dalcin 
7516d9e3a6SLisandro Dalcin #undef __FUNCT__
7616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetUp_HYPRE"
7716d9e3a6SLisandro Dalcin static PetscErrorCode PCSetUp_HYPRE(PC pc)
7816d9e3a6SLisandro Dalcin {
7916d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
8016d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
8116d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
8216d9e3a6SLisandro Dalcin   HYPRE_ParVector    bv,xv;
8316d9e3a6SLisandro Dalcin   PetscInt           bs;
8416d9e3a6SLisandro Dalcin 
8516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
8616d9e3a6SLisandro Dalcin   if (!jac->hypre_type) {
8702a17cd4SBarry Smith     ierr = PCHYPRESetType(pc,"boomeramg");CHKERRQ(ierr);
8816d9e3a6SLisandro Dalcin   }
895f5c5b43SBarry Smith 
905f5c5b43SBarry Smith   if (pc->setupcalled) {
915f5c5b43SBarry Smith     /* always destroy the old matrix and create a new memory;
925f5c5b43SBarry Smith        hope this does not churn the memory too much. The problem
935f5c5b43SBarry Smith        is I do not know if it is possible to put the matrix back to
945f5c5b43SBarry Smith        its initial state so that we can directly copy the values
955f5c5b43SBarry Smith        the second time through. */
96*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_IJMatrixDestroy(jac->ij));
975f5c5b43SBarry Smith     jac->ij = 0;
9816d9e3a6SLisandro Dalcin   }
995f5c5b43SBarry Smith 
10016d9e3a6SLisandro Dalcin   if (!jac->ij) { /* create the matrix the first time through */
10116d9e3a6SLisandro Dalcin     ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr);
10216d9e3a6SLisandro Dalcin   }
10316d9e3a6SLisandro Dalcin   if (!jac->b) { /* create the vectors the first time through */
10416d9e3a6SLisandro Dalcin     Vec x,b;
10516d9e3a6SLisandro Dalcin     ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr);
10616d9e3a6SLisandro Dalcin     ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr);
10716d9e3a6SLisandro Dalcin     ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr);
10816d9e3a6SLisandro Dalcin     ierr = VecDestroy(x);CHKERRQ(ierr);
10916d9e3a6SLisandro Dalcin     ierr = VecDestroy(b);CHKERRQ(ierr);
11016d9e3a6SLisandro Dalcin   }
1115f5c5b43SBarry Smith 
11216d9e3a6SLisandro Dalcin   /* special case for BoomerAMG */
11316d9e3a6SLisandro Dalcin   if (jac->setup == HYPRE_BoomerAMGSetup) {
11416d9e3a6SLisandro Dalcin     ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr);
11516d9e3a6SLisandro Dalcin     if (bs > 1) {
116*a65764d7SBarry Smith       PetscStackCallHypre("",HYPRE_BoomerAMGSetNumFunctions(jac->hsolver,bs));
11716d9e3a6SLisandro Dalcin     }
11816d9e3a6SLisandro Dalcin   };
11916d9e3a6SLisandro Dalcin   ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr);
120*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat));
121*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_IJVectorGetObject(jac->b,(void**)&bv));
122*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_IJVectorGetObject(jac->x,(void**)&xv));
123*a65764d7SBarry Smith   PetscStackCallHypre("HYPRE_SetupXXX",(*jac->setup)(jac->hsolver,hmat,bv,xv));
12416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
12516d9e3a6SLisandro Dalcin }
12616d9e3a6SLisandro Dalcin 
12716d9e3a6SLisandro Dalcin /*
12816d9e3a6SLisandro Dalcin     Replaces the address where the HYPRE vector points to its data with the address of
12916d9e3a6SLisandro Dalcin   PETSc's data. Saves the old address so it can be reset when we are finished with it.
13016d9e3a6SLisandro Dalcin   Allows use to get the data into a HYPRE vector without the cost of memcopies
13116d9e3a6SLisandro Dalcin */
13216d9e3a6SLisandro Dalcin #define HYPREReplacePointer(b,newvalue,savedvalue) {\
13316d9e3a6SLisandro Dalcin    hypre_ParVector *par_vector   = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));\
13416d9e3a6SLisandro Dalcin    hypre_Vector    *local_vector = hypre_ParVectorLocalVector(par_vector);\
13516d9e3a6SLisandro Dalcin    savedvalue         = local_vector->data;\
13616d9e3a6SLisandro Dalcin    local_vector->data = newvalue;}
13716d9e3a6SLisandro Dalcin 
13816d9e3a6SLisandro Dalcin #undef __FUNCT__
13916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApply_HYPRE"
14016d9e3a6SLisandro Dalcin static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x)
14116d9e3a6SLisandro Dalcin {
14216d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
14316d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
14416d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
14516d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
14616d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
14716d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
14816d9e3a6SLisandro Dalcin   int                hierr;
14916d9e3a6SLisandro Dalcin 
15016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
15116d9e3a6SLisandro Dalcin   if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);}
15216d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
15316d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
15416d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
15516d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
15616d9e3a6SLisandro Dalcin 
157*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat));
158*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_IJVectorGetObject(jac->b,(void**)&jbv));
159*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_IJVectorGetObject(jac->x,(void**)&jxv));
16016d9e3a6SLisandro Dalcin   hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv);
16165e19b50SBarry Smith   if (hierr && hierr != HYPRE_ERROR_CONV) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
16216d9e3a6SLisandro Dalcin   if (hierr) hypre__global_error = 0;
16316d9e3a6SLisandro Dalcin 
16416d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
16516d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
16616d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
16716d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
16816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
16916d9e3a6SLisandro Dalcin }
17016d9e3a6SLisandro Dalcin 
17116d9e3a6SLisandro Dalcin #undef __FUNCT__
17216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCDestroy_HYPRE"
17316d9e3a6SLisandro Dalcin static PetscErrorCode PCDestroy_HYPRE(PC pc)
17416d9e3a6SLisandro Dalcin {
17516d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
17616d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
17716d9e3a6SLisandro Dalcin 
17816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
179*a65764d7SBarry Smith   if (jac->ij) { PetscStackCallHypre("",HYPRE_IJMatrixDestroy(jac->ij)); }
180*a65764d7SBarry Smith   if (jac->b)  { PetscStackCallHypre("",HYPRE_IJVectorDestroy(jac->b));  }
181*a65764d7SBarry Smith   if (jac->x)  { PetscStackCallHypre("",HYPRE_IJVectorDestroy(jac->x));  }
182*a65764d7SBarry Smith   if (jac->destroy) { PetscStackCallHypre("HYPRE_DistroyXXX",(*jac->destroy)(jac->hsolver));}
183503cfb0cSBarry Smith   ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr);
18416d9e3a6SLisandro Dalcin   if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);}
18516d9e3a6SLisandro Dalcin   ierr = PetscFree(jac);CHKERRQ(ierr);
18616d9e3a6SLisandro Dalcin 
18716d9e3a6SLisandro Dalcin   ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr);
18816d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);CHKERRQ(ierr);
18916d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);CHKERRQ(ierr);
19016d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
19116d9e3a6SLisandro Dalcin }
19216d9e3a6SLisandro Dalcin 
19316d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
19416d9e3a6SLisandro Dalcin #undef __FUNCT__
19516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut"
19616d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc)
19716d9e3a6SLisandro Dalcin {
19816d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
19916d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
20016d9e3a6SLisandro Dalcin   PetscTruth     flag;
20116d9e3a6SLisandro Dalcin 
20216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
20316d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr);
20416d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr);
20516d9e3a6SLisandro Dalcin   if (flag) {
206*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParCSRPilutSetMaxIter(jac->hsolver,jac->maxiter));
20716d9e3a6SLisandro Dalcin   }
20816d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr);
20916d9e3a6SLisandro Dalcin   if (flag) {
210*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParCSRPilutSetDropTolerance(jac->hsolver,jac->tol));
21116d9e3a6SLisandro Dalcin   }
21216d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr);
21316d9e3a6SLisandro Dalcin   if (flag) {
214*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParCSRPilutSetFactorRowSize(jac->hsolver,jac->factorrowsize));
21516d9e3a6SLisandro Dalcin   }
21616d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
21716d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
21816d9e3a6SLisandro Dalcin }
21916d9e3a6SLisandro Dalcin 
22016d9e3a6SLisandro Dalcin #undef __FUNCT__
22116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Pilut"
22216d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer)
22316d9e3a6SLisandro Dalcin {
22416d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
22516d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
22616d9e3a6SLisandro Dalcin   PetscTruth     iascii;
22716d9e3a6SLisandro Dalcin 
22816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
2292692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
23016d9e3a6SLisandro Dalcin   if (iascii) {
23116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut preconditioning\n");CHKERRQ(ierr);
23216d9e3a6SLisandro Dalcin     if (jac->maxiter != PETSC_DEFAULT) {
23316d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr);
23416d9e3a6SLisandro Dalcin     } else {
23516d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr);
23616d9e3a6SLisandro Dalcin     }
23716d9e3a6SLisandro Dalcin     if (jac->tol != PETSC_DEFAULT) {
23816d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr);
23916d9e3a6SLisandro Dalcin     } else {
24016d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr);
24116d9e3a6SLisandro Dalcin     }
24216d9e3a6SLisandro Dalcin     if (jac->factorrowsize != PETSC_DEFAULT) {
24316d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr);
24416d9e3a6SLisandro Dalcin     } else {
24516d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Pilut: default factor row size \n");CHKERRQ(ierr);
24616d9e3a6SLisandro Dalcin     }
24716d9e3a6SLisandro Dalcin   }
24816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
24916d9e3a6SLisandro Dalcin }
25016d9e3a6SLisandro Dalcin 
25116d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
25216d9e3a6SLisandro Dalcin #undef __FUNCT__
25316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid"
25416d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc)
25516d9e3a6SLisandro Dalcin {
25616d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
25716d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
25816d9e3a6SLisandro Dalcin   PetscTruth     flag;
259390e7148SBarry Smith   char           *args[8],levels[16];
260390e7148SBarry Smith   PetscInt       cnt = 0;
26116d9e3a6SLisandro Dalcin 
26216d9e3a6SLisandro Dalcin   PetscFunctionBegin;
26316d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr);
26416d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr);
26516d9e3a6SLisandro Dalcin   if (flag) {
26665e19b50SBarry Smith     if (jac->levels < 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels);
26716d9e3a6SLisandro Dalcin     sprintf(levels,"%d",jac->levels);
268390e7148SBarry Smith     args[cnt++] = (char*)"-level"; args[cnt++] = levels;
26916d9e3a6SLisandro Dalcin   }
27016d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);CHKERRQ(ierr);
27116d9e3a6SLisandro Dalcin   if (jac->bjilu) {
272390e7148SBarry Smith     args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1";
27316d9e3a6SLisandro Dalcin   }
27416d9e3a6SLisandro Dalcin 
27516d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);CHKERRQ(ierr);
27616d9e3a6SLisandro Dalcin   if (jac->printstatistics) {
277390e7148SBarry Smith     args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1";
278390e7148SBarry Smith     args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1";
27916d9e3a6SLisandro Dalcin   }
28016d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
281390e7148SBarry Smith   if (cnt) {
282*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_EuclidSetParams(jac->hsolver,cnt,args));
283390e7148SBarry Smith   }
28416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
28516d9e3a6SLisandro Dalcin }
28616d9e3a6SLisandro Dalcin 
28716d9e3a6SLisandro Dalcin #undef __FUNCT__
28816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Euclid"
28916d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer)
29016d9e3a6SLisandro Dalcin {
29116d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
29216d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
29316d9e3a6SLisandro Dalcin   PetscTruth     iascii;
29416d9e3a6SLisandro Dalcin 
29516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
2962692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
29716d9e3a6SLisandro Dalcin   if (iascii) {
29816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid preconditioning\n");CHKERRQ(ierr);
29916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr);
30016d9e3a6SLisandro Dalcin     if (jac->bjilu) {
30116d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr);
30216d9e3a6SLisandro Dalcin     }
30316d9e3a6SLisandro Dalcin   }
30416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
30516d9e3a6SLisandro Dalcin }
30616d9e3a6SLisandro Dalcin 
30716d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
30816d9e3a6SLisandro Dalcin 
30916d9e3a6SLisandro Dalcin #undef __FUNCT__
31016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG"
31116d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x)
31216d9e3a6SLisandro Dalcin {
31316d9e3a6SLisandro Dalcin   PC_HYPRE           *jac = (PC_HYPRE*)pc->data;
31416d9e3a6SLisandro Dalcin   PetscErrorCode     ierr;
31516d9e3a6SLisandro Dalcin   HYPRE_ParCSRMatrix hmat;
31616d9e3a6SLisandro Dalcin   PetscScalar        *bv,*xv;
31716d9e3a6SLisandro Dalcin   HYPRE_ParVector    jbv,jxv;
31816d9e3a6SLisandro Dalcin   PetscScalar        *sbv,*sxv;
31916d9e3a6SLisandro Dalcin   int                hierr;
32016d9e3a6SLisandro Dalcin 
32116d9e3a6SLisandro Dalcin   PetscFunctionBegin;
32216d9e3a6SLisandro Dalcin   ierr = VecSet(x,0.0);CHKERRQ(ierr);
32316d9e3a6SLisandro Dalcin   ierr = VecGetArray(b,&bv);CHKERRQ(ierr);
32416d9e3a6SLisandro Dalcin   ierr = VecGetArray(x,&xv);CHKERRQ(ierr);
32516d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,bv,sbv);
32616d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,xv,sxv);
32716d9e3a6SLisandro Dalcin 
328*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat));
329*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_IJVectorGetObject(jac->b,(void**)&jbv));
330*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_IJVectorGetObject(jac->x,(void**)&jxv));
33116d9e3a6SLisandro Dalcin 
33216d9e3a6SLisandro Dalcin   hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv);
33316d9e3a6SLisandro Dalcin   /* error code of 1 in BoomerAMG merely means convergence not achieved */
334e32f2f54SBarry Smith   if (hierr && (hierr != 1)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr);
33516d9e3a6SLisandro Dalcin   if (hierr) hypre__global_error = 0;
33616d9e3a6SLisandro Dalcin 
33716d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->b,sbv,bv);
33816d9e3a6SLisandro Dalcin   HYPREReplacePointer(jac->x,sxv,xv);
33916d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr);
34016d9e3a6SLisandro Dalcin   ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr);
34116d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
34216d9e3a6SLisandro Dalcin }
34316d9e3a6SLisandro Dalcin 
34416d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCycleType[]   = {"","V","W"};
3450f1074feSSatish Balay static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"};
34616d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGMeasureType[] = {"local","global"};
3470f1074feSSatish Balay static const char *HYPREBoomerAMGRelaxType[]   = {"Jacobi","sequential-Gauss-Seidel","","SOR/Jacobi","backward-SOR/Jacobi","","symmetric-SOR/Jacobi",
3480f1074feSSatish Balay                                                   "","","Gaussian-elimination"};
3490f1074feSSatish Balay static const char *HYPREBoomerAMGInterpType[]  = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i",
3500f1074feSSatish Balay                                                   "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"};
35116d9e3a6SLisandro Dalcin #undef __FUNCT__
35216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG"
35316d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc)
35416d9e3a6SLisandro Dalcin {
35516d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
35616d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
35716d9e3a6SLisandro Dalcin   int            n,indx;
35816d9e3a6SLisandro Dalcin   PetscTruth     flg, tmp_truth;
35916d9e3a6SLisandro Dalcin   double         tmpdbl, twodbl[2];
36016d9e3a6SLisandro Dalcin 
36116d9e3a6SLisandro Dalcin   PetscFunctionBegin;
36216d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr);
36316d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr);
36416d9e3a6SLisandro Dalcin   if (flg) {
36516d9e3a6SLisandro Dalcin     jac->cycletype = indx;
366*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype));
36716d9e3a6SLisandro Dalcin   }
36816d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr);
36916d9e3a6SLisandro Dalcin   if (flg) {
37065e19b50SBarry Smith     if (jac->maxlevels < 2) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels);
371*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels));
37216d9e3a6SLisandro Dalcin   }
37316d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr);
37416d9e3a6SLisandro Dalcin   if (flg) {
37565e19b50SBarry Smith     if (jac->maxiter < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter);
376*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter));
37716d9e3a6SLisandro Dalcin   }
3780f1074feSSatish 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);
37916d9e3a6SLisandro Dalcin   if (flg) {
38065e19b50SBarry 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);
381*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol));
38216d9e3a6SLisandro Dalcin   }
38316d9e3a6SLisandro Dalcin 
3840f1074feSSatish Balay   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr);
38516d9e3a6SLisandro Dalcin   if (flg) {
38665e19b50SBarry 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);
387*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor));
38816d9e3a6SLisandro Dalcin   }
38916d9e3a6SLisandro Dalcin 
3900f1074feSSatish 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);
3910f1074feSSatish Balay   if (flg) {
39265e19b50SBarry 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);
393*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax));
3940f1074feSSatish Balay   }
3950f1074feSSatish Balay 
3960f1074feSSatish Balay  ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);CHKERRQ(ierr);
3970f1074feSSatish Balay   if (flg) {
39865e19b50SBarry 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);
3990f1074feSSatish Balay 
400*a65764d7SBarry Smith      PetscStackCallHypre("",HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl));
4010f1074feSSatish Balay   }
4020f1074feSSatish Balay 
4030f1074feSSatish Balay 
4040f1074feSSatish 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);
4050f1074feSSatish Balay   if (flg) {
40665e19b50SBarry 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);
4070f1074feSSatish Balay 
408*a65764d7SBarry Smith      PetscStackCallHypre("",HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths));
4090f1074feSSatish Balay   }
4100f1074feSSatish Balay 
4110f1074feSSatish Balay 
41216d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr);
41316d9e3a6SLisandro Dalcin   if (flg) {
41465e19b50SBarry 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);
415*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold));
41616d9e3a6SLisandro Dalcin   }
41716d9e3a6SLisandro Dalcin   ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr);
41816d9e3a6SLisandro Dalcin   if (flg) {
41965e19b50SBarry 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);
42065e19b50SBarry 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);
421*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum));
42216d9e3a6SLisandro Dalcin   }
42316d9e3a6SLisandro Dalcin 
42416d9e3a6SLisandro Dalcin   /* Grid sweeps */
4250f1074feSSatish 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);
42616d9e3a6SLisandro Dalcin   if (flg) {
427*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetNumSweeps(jac->hsolver,indx));
42816d9e3a6SLisandro Dalcin     /* modify the jac structure so we can view the updated options with PC_View */
42916d9e3a6SLisandro Dalcin     jac->gridsweeps[0] = indx;
4300f1074feSSatish Balay     jac->gridsweeps[1] = indx;
4310f1074feSSatish Balay     /*defaults coarse to 1 */
4320f1074feSSatish Balay     jac->gridsweeps[2] = 1;
43316d9e3a6SLisandro Dalcin   }
4340f1074feSSatish Balay 
4350f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx ,&flg);CHKERRQ(ierr);
43616d9e3a6SLisandro Dalcin   if (flg) {
437*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 1));
4380f1074feSSatish Balay     jac->gridsweeps[0] = indx;
43916d9e3a6SLisandro Dalcin   }
44016d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr);
44116d9e3a6SLisandro Dalcin   if (flg) {
442*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 2));
4430f1074feSSatish Balay     jac->gridsweeps[1] = indx;
44416d9e3a6SLisandro Dalcin   }
4450f1074feSSatish Balay   ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr);
44616d9e3a6SLisandro Dalcin   if (flg) {
447*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 3));
4480f1074feSSatish Balay     jac->gridsweeps[2] = indx;
44916d9e3a6SLisandro Dalcin   }
45016d9e3a6SLisandro Dalcin 
45116d9e3a6SLisandro Dalcin   /* Relax type */
4520f1074feSSatish 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);
45316d9e3a6SLisandro Dalcin   if (flg) {
4540f1074feSSatish Balay     jac->relaxtype[0] = jac->relaxtype[1]  = indx;
455*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetRelaxType(jac->hsolver, indx));
4560f1074feSSatish Balay     /* by default, coarse type set to 9 */
4570f1074feSSatish Balay     jac->relaxtype[2] = 9;
4580f1074feSSatish Balay 
45916d9e3a6SLisandro Dalcin   }
4600f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
46116d9e3a6SLisandro Dalcin   if (flg) {
46216d9e3a6SLisandro Dalcin     jac->relaxtype[0] = indx;
463*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 1));
46416d9e3a6SLisandro Dalcin   }
4650f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr);
46616d9e3a6SLisandro Dalcin   if (flg) {
4670f1074feSSatish Balay     jac->relaxtype[1] = indx;
468*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 2));
46916d9e3a6SLisandro Dalcin   }
47016d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr);
47116d9e3a6SLisandro Dalcin   if (flg) {
4720f1074feSSatish Balay     jac->relaxtype[2] = indx;
473*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 3));
47416d9e3a6SLisandro Dalcin   }
47516d9e3a6SLisandro Dalcin 
47616d9e3a6SLisandro Dalcin   /* Relaxation Weight */
47716d9e3a6SLisandro 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);
47816d9e3a6SLisandro Dalcin   if (flg) {
479*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetRelaxWt(jac->hsolver,tmpdbl));
48016d9e3a6SLisandro Dalcin     jac->relaxweight = tmpdbl;
48116d9e3a6SLisandro Dalcin   }
48216d9e3a6SLisandro Dalcin 
48316d9e3a6SLisandro Dalcin   n=2;
48416d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
48516d9e3a6SLisandro 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);
48616d9e3a6SLisandro Dalcin   if (flg) {
48716d9e3a6SLisandro Dalcin     if (n == 2) {
48816d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
489*a65764d7SBarry Smith       PetscStackCallHypre("",HYPRE_BoomerAMGSetLevelRelaxWt(jac->hsolver,twodbl[0],indx));
49016d9e3a6SLisandro Dalcin     } else {
49165e19b50SBarry 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);
49216d9e3a6SLisandro Dalcin     }
49316d9e3a6SLisandro Dalcin   }
49416d9e3a6SLisandro Dalcin 
49516d9e3a6SLisandro Dalcin   /* Outer relaxation Weight */
49616d9e3a6SLisandro 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);
49716d9e3a6SLisandro Dalcin   if (flg) {
498*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetOuterWt( jac->hsolver, tmpdbl));
49916d9e3a6SLisandro Dalcin     jac->outerrelaxweight = tmpdbl;
50016d9e3a6SLisandro Dalcin   }
50116d9e3a6SLisandro Dalcin 
50216d9e3a6SLisandro Dalcin   n=2;
50316d9e3a6SLisandro Dalcin   twodbl[0] = twodbl[1] = 1.0;
50416d9e3a6SLisandro 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);
50516d9e3a6SLisandro Dalcin   if (flg) {
50616d9e3a6SLisandro Dalcin     if (n == 2) {
50716d9e3a6SLisandro Dalcin       indx =  (int)PetscAbsReal(twodbl[1]);
508*a65764d7SBarry Smith       PetscStackCallHypre("",HYPRE_BoomerAMGSetLevelOuterWt( jac->hsolver, twodbl[0], indx));
50916d9e3a6SLisandro Dalcin     } else {
51065e19b50SBarry 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);
51116d9e3a6SLisandro Dalcin     }
51216d9e3a6SLisandro Dalcin   }
51316d9e3a6SLisandro Dalcin 
51416d9e3a6SLisandro Dalcin   /* the Relax Order */
51516d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
51616d9e3a6SLisandro Dalcin 
51716d9e3a6SLisandro Dalcin   if (flg) {
51816d9e3a6SLisandro Dalcin     jac->relaxorder = 0;
519*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder));
52016d9e3a6SLisandro Dalcin   }
52116d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr);
52216d9e3a6SLisandro Dalcin   if (flg) {
52316d9e3a6SLisandro Dalcin     jac->measuretype = indx;
524*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype));
52516d9e3a6SLisandro Dalcin   }
5260f1074feSSatish Balay   /* update list length 3/07 */
5270f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,11,HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr);
52816d9e3a6SLisandro Dalcin   if (flg) {
52916d9e3a6SLisandro Dalcin     jac->coarsentype = indx;
530*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype));
53116d9e3a6SLisandro Dalcin   }
5320f1074feSSatish Balay 
5330f1074feSSatish Balay   /* new 3/07 */
5340f1074feSSatish Balay   ierr = PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,14,HYPREBoomerAMGInterpType[0],&indx,&flg);CHKERRQ(ierr);
5350f1074feSSatish Balay   if (flg) {
5360f1074feSSatish Balay     jac->interptype = indx;
537*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype));
5380f1074feSSatish Balay   }
5390f1074feSSatish Balay 
54090d69ab7SBarry Smith   flg  = PETSC_FALSE;
54190d69ab7SBarry Smith   ierr = PetscOptionsTruth("-pc_hypre_boomeramg_print_statistics","Print statistics","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
54216d9e3a6SLisandro Dalcin   if (flg) {
54316d9e3a6SLisandro Dalcin     int level=3;
54416d9e3a6SLisandro Dalcin     jac->printstatistics = PETSC_TRUE;
54516d9e3a6SLisandro Dalcin     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
546*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetPrintLevel(jac->hsolver,level));
5472ae77aedSBarry Smith   }
5482ae77aedSBarry Smith 
54990d69ab7SBarry Smith   flg  = PETSC_FALSE;
55090d69ab7SBarry Smith   ierr = PetscOptionsTruth("-pc_hypre_boomeramg_print_debug","Print debug information","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
5512ae77aedSBarry Smith   if (flg) {
5522ae77aedSBarry Smith     int level=3;
5532ae77aedSBarry Smith     jac->printstatistics = PETSC_TRUE;
5542ae77aedSBarry Smith     ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,PETSC_NULL);CHKERRQ(ierr);
555*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetDebugFlag(jac->hsolver,level));
55616d9e3a6SLisandro Dalcin   }
5578f87f92bSBarry Smith 
5588f87f92bSBarry Smith   ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5598f87f92bSBarry Smith   if (flg && tmp_truth) {
5608f87f92bSBarry Smith     jac->nodal_coarsen = 1;
561*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetNodal(jac->hsolver,1));
5628f87f92bSBarry Smith   }
5638f87f92bSBarry Smith 
5648f87f92bSBarry Smith   ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr);
5658f87f92bSBarry Smith   if (flg && tmp_truth) {
5668f87f92bSBarry Smith     PetscInt tmp_int;
5678f87f92bSBarry Smith     ierr = PetscOptionsInt( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);CHKERRQ(ierr);
5688f87f92bSBarry Smith     if (flg) jac->nodal_relax_levels = tmp_int;
569*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetSmoothType(jac->hsolver,6));
570*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetDomainType(jac->hsolver,1));
571*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetOverlap(jac->hsolver,0));
572*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetSmoothNumLevels(jac->hsolver,jac->nodal_relax_levels));
5738f87f92bSBarry Smith   }
5748f87f92bSBarry Smith 
57516d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
57616d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
57716d9e3a6SLisandro Dalcin }
57816d9e3a6SLisandro Dalcin 
57916d9e3a6SLisandro Dalcin #undef __FUNCT__
58016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG"
5817319c654SBarry Smith static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
58216d9e3a6SLisandro Dalcin {
58316d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
58416d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
5854d0a8057SBarry Smith   int            oits;
58616d9e3a6SLisandro Dalcin 
58716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
588*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_BoomerAMGSetMaxIter(jac->hsolver,its*jac->maxiter));
589*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_BoomerAMGSetTol(jac->hsolver,rtol));
59016d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_TRUE;
59116d9e3a6SLisandro Dalcin   ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr);
59216d9e3a6SLisandro Dalcin   jac->applyrichardson = PETSC_FALSE;
593*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_BoomerAMGGetNumIterations(jac->hsolver,&oits));
5944d0a8057SBarry Smith   *outits = oits;
5954d0a8057SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
5964d0a8057SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
597*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol));
598*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter));
59916d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
60016d9e3a6SLisandro Dalcin }
60116d9e3a6SLisandro Dalcin 
60216d9e3a6SLisandro Dalcin 
60316d9e3a6SLisandro Dalcin #undef __FUNCT__
60416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_BoomerAMG"
60516d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer)
60616d9e3a6SLisandro Dalcin {
60716d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
60816d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
60916d9e3a6SLisandro Dalcin   PetscTruth     iascii;
61016d9e3a6SLisandro Dalcin 
61116d9e3a6SLisandro Dalcin   PetscFunctionBegin;
6122692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
61316d9e3a6SLisandro Dalcin   if (iascii) {
61416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr);
61516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr);
61616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr);
61716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr);
61816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr);
61916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr);
6200f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation truncation factor %G\n",jac->truncfactor);CHKERRQ(ierr);
6210f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);CHKERRQ(ierr);
6220f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);CHKERRQ(ierr);
6230f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);CHKERRQ(ierr);
6240f1074feSSatish Balay 
62516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr);
62616d9e3a6SLisandro Dalcin 
6270f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps down         %d\n",jac->gridsweeps[0]);CHKERRQ(ierr);
6280f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps up           %d\n",jac->gridsweeps[1]);CHKERRQ(ierr);
6290f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Sweeps on coarse    %d\n",jac->gridsweeps[2]);CHKERRQ(ierr);
63016d9e3a6SLisandro Dalcin 
6310f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax down          %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr);
6320f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax up            %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr);
6330f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax on coarse     %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr);
63416d9e3a6SLisandro Dalcin 
63516d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Relax weight  (all)      %G\n",jac->relaxweight);CHKERRQ(ierr);
63616d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr);
63716d9e3a6SLisandro Dalcin 
63816d9e3a6SLisandro Dalcin     if (jac->relaxorder) {
63916d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr);
64016d9e3a6SLisandro Dalcin     } else {
64116d9e3a6SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr);
64216d9e3a6SLisandro Dalcin     }
64316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Measure type        %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr);
64416d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Coarsen type        %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr);
6450f1074feSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE BoomerAMG: Interpolation type  %s\n",HYPREBoomerAMGInterpType[jac->interptype]);CHKERRQ(ierr);
6468f87f92bSBarry Smith     if (jac->nodal_coarsen) {
6478f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal())\n");CHKERRQ(ierr);
6488f87f92bSBarry Smith     }
6498f87f92bSBarry Smith     if (jac->nodal_relax) {
6508f87f92bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);CHKERRQ(ierr);
6518f87f92bSBarry Smith     }
65216d9e3a6SLisandro Dalcin   }
65316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
65416d9e3a6SLisandro Dalcin }
65516d9e3a6SLisandro Dalcin 
65616d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/
65716d9e3a6SLisandro Dalcin #undef __FUNCT__
65816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails"
65916d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc)
66016d9e3a6SLisandro Dalcin {
66116d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
66216d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
66316d9e3a6SLisandro Dalcin   int            indx;
66416d9e3a6SLisandro Dalcin   PetscTruth     flag;
66516d9e3a6SLisandro Dalcin   const char     *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"};
66616d9e3a6SLisandro Dalcin 
66716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
66816d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr);
66916d9e3a6SLisandro Dalcin   ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr);
67016d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr);
67116d9e3a6SLisandro Dalcin   if (flag) {
672*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels));
67316d9e3a6SLisandro Dalcin   }
67416d9e3a6SLisandro Dalcin 
67516d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr);
67616d9e3a6SLisandro Dalcin   if (flag) {
677*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter));
67816d9e3a6SLisandro Dalcin   }
67916d9e3a6SLisandro Dalcin 
68016d9e3a6SLisandro Dalcin   ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr);
68116d9e3a6SLisandro Dalcin   if (flag) {
682*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal));
68316d9e3a6SLisandro Dalcin   }
68416d9e3a6SLisandro Dalcin 
68516d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,&flag);CHKERRQ(ierr);
68616d9e3a6SLisandro Dalcin   if (flag) {
687*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging));
68816d9e3a6SLisandro Dalcin   }
68916d9e3a6SLisandro Dalcin 
69016d9e3a6SLisandro Dalcin   ierr = PetscOptionsTruth("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,&flag);CHKERRQ(ierr);
69116d9e3a6SLisandro Dalcin   if (flag) {
692*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse));
69316d9e3a6SLisandro Dalcin   }
69416d9e3a6SLisandro Dalcin 
69516d9e3a6SLisandro Dalcin   ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);CHKERRQ(ierr);
69616d9e3a6SLisandro Dalcin   if (flag) {
69716d9e3a6SLisandro Dalcin     jac->symt = indx;
698*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt));
69916d9e3a6SLisandro Dalcin   }
70016d9e3a6SLisandro Dalcin 
70116d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
70216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
70316d9e3a6SLisandro Dalcin }
70416d9e3a6SLisandro Dalcin 
70516d9e3a6SLisandro Dalcin #undef __FUNCT__
70616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_ParaSails"
70716d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer)
70816d9e3a6SLisandro Dalcin {
70916d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
71016d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
71116d9e3a6SLisandro Dalcin   PetscTruth     iascii;
71216d9e3a6SLisandro Dalcin   const char     *symt = 0;;
71316d9e3a6SLisandro Dalcin 
71416d9e3a6SLisandro Dalcin   PetscFunctionBegin;
7152692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
71616d9e3a6SLisandro Dalcin   if (iascii) {
71716d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails preconditioning\n");CHKERRQ(ierr);
71816d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr);
71916d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr);
72016d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr);
72116d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr);
72216d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: reuse nonzero structure %s\n",PetscTruths[jac->ruse]);CHKERRQ(ierr);
72316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: print info to screen %s\n",PetscTruths[jac->logging]);CHKERRQ(ierr);
72416d9e3a6SLisandro Dalcin     if (!jac->symt) {
72516d9e3a6SLisandro Dalcin       symt = "nonsymmetric matrix and preconditioner";
72616d9e3a6SLisandro Dalcin     } else if (jac->symt == 1) {
72716d9e3a6SLisandro Dalcin       symt = "SPD matrix and preconditioner";
72816d9e3a6SLisandro Dalcin     } else if (jac->symt == 2) {
72916d9e3a6SLisandro Dalcin       symt = "nonsymmetric matrix but SPD preconditioner";
73016d9e3a6SLisandro Dalcin     } else {
73165e19b50SBarry Smith       SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt);
73216d9e3a6SLisandro Dalcin     }
73316d9e3a6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr);
73416d9e3a6SLisandro Dalcin   }
73516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
73616d9e3a6SLisandro Dalcin }
73716d9e3a6SLisandro Dalcin /* ---------------------------------------------------------------------------------*/
73816d9e3a6SLisandro Dalcin 
73916d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
74016d9e3a6SLisandro Dalcin #undef __FUNCT__
74116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType_HYPRE"
74216d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType_HYPRE(PC pc,const char *name[])
74316d9e3a6SLisandro Dalcin {
74416d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
74516d9e3a6SLisandro Dalcin 
74616d9e3a6SLisandro Dalcin   PetscFunctionBegin;
74716d9e3a6SLisandro Dalcin   *name = jac->hypre_type;
74816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
74916d9e3a6SLisandro Dalcin }
75016d9e3a6SLisandro Dalcin EXTERN_C_END
75116d9e3a6SLisandro Dalcin 
75216d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
75316d9e3a6SLisandro Dalcin #undef __FUNCT__
75416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType_HYPRE"
75516d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType_HYPRE(PC pc,const char name[])
75616d9e3a6SLisandro Dalcin {
75716d9e3a6SLisandro Dalcin   PC_HYPRE       *jac = (PC_HYPRE*)pc->data;
75816d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
75916d9e3a6SLisandro Dalcin   PetscTruth     flag;
76016d9e3a6SLisandro Dalcin 
76116d9e3a6SLisandro Dalcin   PetscFunctionBegin;
76216d9e3a6SLisandro Dalcin   if (jac->hypre_type) {
76316d9e3a6SLisandro Dalcin     ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr);
764e7e72b3dSBarry Smith     if (!flag) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set");
76516d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
76616d9e3a6SLisandro Dalcin   } else {
76716d9e3a6SLisandro Dalcin     ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr);
76816d9e3a6SLisandro Dalcin   }
76916d9e3a6SLisandro Dalcin 
77016d9e3a6SLisandro Dalcin   jac->maxiter            = PETSC_DEFAULT;
77116d9e3a6SLisandro Dalcin   jac->tol                = PETSC_DEFAULT;
77216d9e3a6SLisandro Dalcin   jac->printstatistics    = PetscLogPrintInfo;
77316d9e3a6SLisandro Dalcin 
77416d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr);
77516d9e3a6SLisandro Dalcin   if (flag) {
776*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParCSRPilutCreate(jac->comm_hypre,&jac->hsolver));
77716d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut;
77816d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Pilut;
77916d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParCSRPilutDestroy;
78016d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParCSRPilutSetup;
78116d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParCSRPilutSolve;
78216d9e3a6SLisandro Dalcin     jac->factorrowsize      = PETSC_DEFAULT;
78316d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
78416d9e3a6SLisandro Dalcin   }
78516d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr);
78616d9e3a6SLisandro Dalcin   if (flag) {
787*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParaSailsCreate(jac->comm_hypre,&jac->hsolver));
78816d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails;
78916d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_ParaSails;
79016d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_ParaSailsDestroy;
79116d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_ParaSailsSetup;
79216d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_ParaSailsSolve;
79316d9e3a6SLisandro Dalcin     /* initialize */
79416d9e3a6SLisandro Dalcin     jac->nlevels     = 1;
79516d9e3a6SLisandro Dalcin     jac->threshhold  = .1;
79616d9e3a6SLisandro Dalcin     jac->filter      = .1;
79716d9e3a6SLisandro Dalcin     jac->loadbal     = 0;
79816d9e3a6SLisandro Dalcin     if (PetscLogPrintInfo) {
79916d9e3a6SLisandro Dalcin       jac->logging     = (int) PETSC_TRUE;
80016d9e3a6SLisandro Dalcin     } else {
80116d9e3a6SLisandro Dalcin       jac->logging     = (int) PETSC_FALSE;
80216d9e3a6SLisandro Dalcin     }
80316d9e3a6SLisandro Dalcin     jac->ruse = (int) PETSC_FALSE;
80416d9e3a6SLisandro Dalcin     jac->symt   = 0;
805*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels));
806*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter));
807*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal));
808*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging));
809*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse));
810*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt));
81116d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
81216d9e3a6SLisandro Dalcin   }
81316d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr);
81416d9e3a6SLisandro Dalcin   if (flag) {
81516d9e3a6SLisandro Dalcin     ierr                    = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver);
81616d9e3a6SLisandro Dalcin     pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid;
81716d9e3a6SLisandro Dalcin     pc->ops->view           = PCView_HYPRE_Euclid;
81816d9e3a6SLisandro Dalcin     jac->destroy            = HYPRE_EuclidDestroy;
81916d9e3a6SLisandro Dalcin     jac->setup              = HYPRE_EuclidSetup;
82016d9e3a6SLisandro Dalcin     jac->solve              = HYPRE_EuclidSolve;
82116d9e3a6SLisandro Dalcin     /* initialization */
82216d9e3a6SLisandro Dalcin     jac->bjilu              = PETSC_FALSE;
82316d9e3a6SLisandro Dalcin     jac->levels             = 1;
82416d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
82516d9e3a6SLisandro Dalcin   }
82616d9e3a6SLisandro Dalcin   ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr);
82716d9e3a6SLisandro Dalcin   if (flag) {
82816d9e3a6SLisandro Dalcin     ierr                     = HYPRE_BoomerAMGCreate(&jac->hsolver);
82916d9e3a6SLisandro Dalcin     pc->ops->setfromoptions  = PCSetFromOptions_HYPRE_BoomerAMG;
83016d9e3a6SLisandro Dalcin     pc->ops->view            = PCView_HYPRE_BoomerAMG;
83116d9e3a6SLisandro Dalcin     pc->ops->applytranspose  = PCApplyTranspose_HYPRE_BoomerAMG;
83216d9e3a6SLisandro Dalcin     pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG;
83316d9e3a6SLisandro Dalcin     jac->destroy             = HYPRE_BoomerAMGDestroy;
83416d9e3a6SLisandro Dalcin     jac->setup               = HYPRE_BoomerAMGSetup;
83516d9e3a6SLisandro Dalcin     jac->solve               = HYPRE_BoomerAMGSolve;
83616d9e3a6SLisandro Dalcin     jac->applyrichardson     = PETSC_FALSE;
83716d9e3a6SLisandro Dalcin     /* these defaults match the hypre defaults */
83816d9e3a6SLisandro Dalcin     jac->cycletype        = 1;
83916d9e3a6SLisandro Dalcin     jac->maxlevels        = 25;
84016d9e3a6SLisandro Dalcin     jac->maxiter          = 1;
8418f87f92bSBarry Smith     jac->tol              = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */
84216d9e3a6SLisandro Dalcin     jac->truncfactor      = 0.0;
84316d9e3a6SLisandro Dalcin     jac->strongthreshold  = .25;
84416d9e3a6SLisandro Dalcin     jac->maxrowsum        = .9;
84516d9e3a6SLisandro Dalcin     jac->coarsentype      = 6;
84616d9e3a6SLisandro Dalcin     jac->measuretype      = 0;
8470f1074feSSatish Balay     jac->gridsweeps[0]    = jac->gridsweeps[1] = jac->gridsweeps[2] = 1;
8488f87f92bSBarry Smith     jac->relaxtype[0]     = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */
8490f1074feSSatish Balay     jac->relaxtype[2]     = 9; /*G.E. */
85016d9e3a6SLisandro Dalcin     jac->relaxweight      = 1.0;
85116d9e3a6SLisandro Dalcin     jac->outerrelaxweight = 1.0;
85216d9e3a6SLisandro Dalcin     jac->relaxorder       = 1;
8530f1074feSSatish Balay     jac->interptype       = 0;
8540f1074feSSatish Balay     jac->agg_nl           = 0;
8550f1074feSSatish Balay     jac->pmax             = 0;
8560f1074feSSatish Balay     jac->truncfactor      = 0.0;
8570f1074feSSatish Balay     jac->agg_num_paths    = 1;
8588f87f92bSBarry Smith 
8598f87f92bSBarry Smith     jac->nodal_coarsen    = 0;
8608f87f92bSBarry Smith     jac->nodal_relax      = PETSC_FALSE;
8618f87f92bSBarry Smith     jac->nodal_relax_levels = 1;
862*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype));
863*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels));
864*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter));
865*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol));
866*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor));
867*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold));
868*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum));
869*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype));
870*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype));
871*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder));
872*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype));
873*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl));
874*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax));
875*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths));
876*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetRelaxType(jac->hsolver, jac->relaxtype[0]));  /*defaults coarse to 9*/
877*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_BoomerAMGSetNumSweeps(jac->hsolver, jac->gridsweeps[0])); /*defaults coarse to 1 */
87816d9e3a6SLisandro Dalcin     PetscFunctionReturn(0);
87916d9e3a6SLisandro Dalcin   }
880503cfb0cSBarry Smith   ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr);
88116d9e3a6SLisandro Dalcin   jac->hypre_type = PETSC_NULL;
88265e19b50SBarry Smith   SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name);
88316d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
88416d9e3a6SLisandro Dalcin }
88516d9e3a6SLisandro Dalcin EXTERN_C_END
88616d9e3a6SLisandro Dalcin 
88716d9e3a6SLisandro Dalcin /*
88816d9e3a6SLisandro Dalcin     It only gets here if the HYPRE type has not been set before the call to
88916d9e3a6SLisandro Dalcin    ...SetFromOptions() which actually is most of the time
89016d9e3a6SLisandro Dalcin */
89116d9e3a6SLisandro Dalcin #undef __FUNCT__
89216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE"
89316d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE(PC pc)
89416d9e3a6SLisandro Dalcin {
89516d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
89616d9e3a6SLisandro Dalcin   int            indx;
89716d9e3a6SLisandro Dalcin   const char     *type[] = {"pilut","parasails","boomeramg","euclid"};
89816d9e3a6SLisandro Dalcin   PetscTruth     flg;
89916d9e3a6SLisandro Dalcin 
90016d9e3a6SLisandro Dalcin   PetscFunctionBegin;
90116d9e3a6SLisandro Dalcin   ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr);
90202a17cd4SBarry Smith   ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"boomeramg",&indx,&flg);CHKERRQ(ierr);
90316d9e3a6SLisandro Dalcin   if (flg) {
90416d9e3a6SLisandro Dalcin     ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr);
90502a17cd4SBarry Smith   } else {
90602a17cd4SBarry Smith     ierr = PCHYPRESetType_HYPRE(pc,"boomeramg");CHKERRQ(ierr);
90716d9e3a6SLisandro Dalcin   }
90816d9e3a6SLisandro Dalcin   if (pc->ops->setfromoptions) {
90916d9e3a6SLisandro Dalcin     ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr);
91016d9e3a6SLisandro Dalcin   }
91116d9e3a6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
91216d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
91316d9e3a6SLisandro Dalcin }
91416d9e3a6SLisandro Dalcin 
91516d9e3a6SLisandro Dalcin #undef __FUNCT__
91616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType"
91716d9e3a6SLisandro Dalcin /*@C
91816d9e3a6SLisandro Dalcin      PCHYPRESetType - Sets which hypre preconditioner you wish to use
91916d9e3a6SLisandro Dalcin 
92016d9e3a6SLisandro Dalcin    Input Parameters:
92116d9e3a6SLisandro Dalcin +     pc - the preconditioner context
92216d9e3a6SLisandro Dalcin -     name - either  pilut, parasails, boomeramg, euclid
92316d9e3a6SLisandro Dalcin 
92416d9e3a6SLisandro Dalcin    Options Database Keys:
92516d9e3a6SLisandro Dalcin    -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
92616d9e3a6SLisandro Dalcin 
92716d9e3a6SLisandro Dalcin    Level: intermediate
92816d9e3a6SLisandro Dalcin 
92916d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
93016d9e3a6SLisandro Dalcin            PCHYPRE
93116d9e3a6SLisandro Dalcin 
93216d9e3a6SLisandro Dalcin @*/
93316d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC pc,const char name[])
93416d9e3a6SLisandro Dalcin {
93516d9e3a6SLisandro Dalcin   PetscErrorCode ierr,(*f)(PC,const char[]);
93616d9e3a6SLisandro Dalcin 
93716d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9380700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
93916d9e3a6SLisandro Dalcin   PetscValidCharPointer(name,2);
94016d9e3a6SLisandro Dalcin   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPRESetType_C",(void (**)(void))&f);CHKERRQ(ierr);
94116d9e3a6SLisandro Dalcin   if (f) {
94216d9e3a6SLisandro Dalcin     ierr = (*f)(pc,name);CHKERRQ(ierr);
94316d9e3a6SLisandro Dalcin   }
94416d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
94516d9e3a6SLisandro Dalcin }
94616d9e3a6SLisandro Dalcin 
94716d9e3a6SLisandro Dalcin #undef __FUNCT__
94816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType"
94916d9e3a6SLisandro Dalcin /*@C
95016d9e3a6SLisandro Dalcin      PCHYPREGetType - Gets which hypre preconditioner you are using
95116d9e3a6SLisandro Dalcin 
95216d9e3a6SLisandro Dalcin    Input Parameter:
95316d9e3a6SLisandro Dalcin .     pc - the preconditioner context
95416d9e3a6SLisandro Dalcin 
95516d9e3a6SLisandro Dalcin    Output Parameter:
95616d9e3a6SLisandro Dalcin .     name - either  pilut, parasails, boomeramg, euclid
95716d9e3a6SLisandro Dalcin 
95816d9e3a6SLisandro Dalcin    Level: intermediate
95916d9e3a6SLisandro Dalcin 
96016d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC,
96116d9e3a6SLisandro Dalcin            PCHYPRE
96216d9e3a6SLisandro Dalcin 
96316d9e3a6SLisandro Dalcin @*/
96416d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC pc,const char *name[])
96516d9e3a6SLisandro Dalcin {
96616d9e3a6SLisandro Dalcin   PetscErrorCode ierr,(*f)(PC,const char*[]);
96716d9e3a6SLisandro Dalcin 
96816d9e3a6SLisandro Dalcin   PetscFunctionBegin;
9690700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
97016d9e3a6SLisandro Dalcin   PetscValidPointer(name,2);
97116d9e3a6SLisandro Dalcin   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPREGetType_C",(void (**)(void))&f);CHKERRQ(ierr);
97216d9e3a6SLisandro Dalcin   if (f) {
97316d9e3a6SLisandro Dalcin     ierr = (*f)(pc,name);CHKERRQ(ierr);
97416d9e3a6SLisandro Dalcin   }
97516d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
97616d9e3a6SLisandro Dalcin }
97716d9e3a6SLisandro Dalcin 
97816d9e3a6SLisandro Dalcin /*MC
97916d9e3a6SLisandro Dalcin      PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre
98016d9e3a6SLisandro Dalcin 
98116d9e3a6SLisandro Dalcin    Options Database Keys:
98216d9e3a6SLisandro Dalcin +   -pc_hypre_type - One of pilut, parasails, boomeramg, euclid
98316d9e3a6SLisandro Dalcin -   Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX
98416d9e3a6SLisandro Dalcin           preconditioner
98516d9e3a6SLisandro Dalcin 
98616d9e3a6SLisandro Dalcin    Level: intermediate
98716d9e3a6SLisandro Dalcin 
98816d9e3a6SLisandro Dalcin    Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()),
98916d9e3a6SLisandro Dalcin           the many hypre options can ONLY be set via the options database (e.g. the command line
99016d9e3a6SLisandro Dalcin           or with PetscOptionsSetValue(), there are no functions to set them)
99116d9e3a6SLisandro Dalcin 
99216d9e3a6SLisandro Dalcin           The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations
9930f1074feSSatish Balay           (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if
9940f1074feSSatish Balay           -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner
9950f1074feSSatish Balay           (-pc_hypre_boomeramg_rtol should be set to 0.0 - the default - to strictly use a fixed number of
9968f87f92bSBarry Smith           iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations
9970f1074feSSatish Balay           and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10
9980f1074feSSatish Balay           then AT MOST twenty V-cycles of boomeramg will be called.
99916d9e3a6SLisandro Dalcin 
10000f1074feSSatish Balay            Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation
10010f1074feSSatish Balay            (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry.
10020f1074feSSatish Balay            Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi.
100316d9e3a6SLisandro Dalcin           If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly
100416d9e3a6SLisandro Dalcin           and use -ksp_max_it to control the number of V-cycles.
100516d9e3a6SLisandro Dalcin           (see the PETSc FAQ.html at the PETSc website under the Documentation tab).
100616d9e3a6SLisandro Dalcin 
100716d9e3a6SLisandro Dalcin 	  2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option
100816d9e3a6SLisandro Dalcin 	  -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L.
100916d9e3a6SLisandro Dalcin 
10109e5bc791SBarry Smith           See PCPFMG for access to the hypre Struct PFMG solver
10119e5bc791SBarry Smith 
101216d9e3a6SLisandro Dalcin .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
10139e5bc791SBarry Smith            PCHYPRESetType(), PCPFMG
101416d9e3a6SLisandro Dalcin 
101516d9e3a6SLisandro Dalcin M*/
101616d9e3a6SLisandro Dalcin 
101716d9e3a6SLisandro Dalcin EXTERN_C_BEGIN
101816d9e3a6SLisandro Dalcin #undef __FUNCT__
101916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCCreate_HYPRE"
102016d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_HYPRE(PC pc)
102116d9e3a6SLisandro Dalcin {
102216d9e3a6SLisandro Dalcin   PC_HYPRE       *jac;
102316d9e3a6SLisandro Dalcin   PetscErrorCode ierr;
102416d9e3a6SLisandro Dalcin 
102516d9e3a6SLisandro Dalcin   PetscFunctionBegin;
102638f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_HYPRE,&jac);CHKERRQ(ierr);
102716d9e3a6SLisandro Dalcin   pc->data                 = jac;
102816d9e3a6SLisandro Dalcin   pc->ops->destroy         = PCDestroy_HYPRE;
102916d9e3a6SLisandro Dalcin   pc->ops->setfromoptions  = PCSetFromOptions_HYPRE;
103016d9e3a6SLisandro Dalcin   pc->ops->setup           = PCSetUp_HYPRE;
103116d9e3a6SLisandro Dalcin   pc->ops->apply           = PCApply_HYPRE;
103216d9e3a6SLisandro Dalcin   jac->comm_hypre          = MPI_COMM_NULL;
103316d9e3a6SLisandro Dalcin   jac->hypre_type          = PETSC_NULL;
103416d9e3a6SLisandro Dalcin   /* duplicate communicator for hypre */
10357adad957SLisandro Dalcin   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(jac->comm_hypre));CHKERRQ(ierr);
103616d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr);
103716d9e3a6SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr);
103816d9e3a6SLisandro Dalcin   PetscFunctionReturn(0);
103916d9e3a6SLisandro Dalcin }
104016d9e3a6SLisandro Dalcin EXTERN_C_END
1041ebc551c0SBarry Smith 
1042f91d8e95SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------*/
1043f91d8e95SBarry Smith 
1044b862ddfaSBarry Smith /* this include is needed ONLY to allow access to the private data inside the Mat object specific to hypre */
104568326731SBarry Smith #include "private/matimpl.h"
1046ebc551c0SBarry Smith 
1047ebc551c0SBarry Smith typedef struct {
104868326731SBarry Smith   MPI_Comm            hcomm;       /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */
1049f91d8e95SBarry Smith   HYPRE_StructSolver  hsolver;
10509e5bc791SBarry Smith 
10519e5bc791SBarry Smith   /* keep copy of PFMG options used so may view them */
10529e5bc791SBarry Smith   int                 its;
10539e5bc791SBarry Smith   double              tol;
10549e5bc791SBarry Smith   int                 relax_type;
10559e5bc791SBarry Smith   int                 rap_type;
10569e5bc791SBarry Smith   int                 num_pre_relax,num_post_relax;
10573b46a515SGlenn Hammond   int                 max_levels;
1058ebc551c0SBarry Smith } PC_PFMG;
1059ebc551c0SBarry Smith 
1060ebc551c0SBarry Smith #undef __FUNCT__
1061ebc551c0SBarry Smith #define __FUNCT__ "PCDestroy_PFMG"
1062ebc551c0SBarry Smith PetscErrorCode PCDestroy_PFMG(PC pc)
1063ebc551c0SBarry Smith {
1064ebc551c0SBarry Smith   PetscErrorCode ierr;
1065f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1066ebc551c0SBarry Smith 
1067ebc551c0SBarry Smith   PetscFunctionBegin;
1068*a65764d7SBarry Smith   if (ex->hsolver) {PetscStackCallHypre("",HYPRE_StructPFMGDestroy(ex->hsolver));}
1069f91d8e95SBarry Smith   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1070f91d8e95SBarry Smith   ierr = PetscFree(ex);CHKERRQ(ierr);
1071ebc551c0SBarry Smith   PetscFunctionReturn(0);
1072ebc551c0SBarry Smith }
1073ebc551c0SBarry Smith 
10749e5bc791SBarry Smith static const char *PFMGRelaxType[]   = {"Jacobi","Weighted-Jacobi","symmetric-Red/Black-Gauss-Seidel","Red/Black-Gauss-Seidel"};
10759e5bc791SBarry Smith static const char *PFMGRAPType[]   = {"Galerkin","non-Galerkin"};
10769e5bc791SBarry Smith 
1077ebc551c0SBarry Smith #undef __FUNCT__
1078ebc551c0SBarry Smith #define __FUNCT__ "PCView_PFMG"
1079ebc551c0SBarry Smith PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer)
1080ebc551c0SBarry Smith {
1081ebc551c0SBarry Smith   PetscErrorCode ierr;
1082ebc551c0SBarry Smith   PetscTruth     iascii;
1083f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1084ebc551c0SBarry Smith 
1085ebc551c0SBarry Smith   PetscFunctionBegin;
10862692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
10879e5bc791SBarry Smith   if (iascii) {
10889e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG preconditioning\n");CHKERRQ(ierr);
10899e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
10909e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
10919e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
10929e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: RAP type %s\n",PFMGRAPType[ex->rap_type]);CHKERRQ(ierr);
10939e5bc791SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
10943b46a515SGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE PFMG: max levels %d\n",ex->max_levels);CHKERRQ(ierr);
10959e5bc791SBarry Smith   }
1096ebc551c0SBarry Smith   PetscFunctionReturn(0);
1097ebc551c0SBarry Smith }
1098ebc551c0SBarry Smith 
10999e5bc791SBarry Smith 
1100ebc551c0SBarry Smith #undef __FUNCT__
1101ebc551c0SBarry Smith #define __FUNCT__ "PCSetFromOptions_PFMG"
1102ebc551c0SBarry Smith PetscErrorCode PCSetFromOptions_PFMG(PC pc)
1103ebc551c0SBarry Smith {
1104ebc551c0SBarry Smith   PetscErrorCode ierr;
1105f91d8e95SBarry Smith   PC_PFMG        *ex = (PC_PFMG*) pc->data;
1106886061d4SBarry Smith   PetscTruth     flg = PETSC_FALSE;
1107ebc551c0SBarry Smith 
1108ebc551c0SBarry Smith   PetscFunctionBegin;
1109ebc551c0SBarry Smith   ierr = PetscOptionsHead("PFMG options");CHKERRQ(ierr);
11109e5bc791SBarry Smith   ierr = PetscOptionsTruth("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
111168326731SBarry Smith   if (flg) {
1112a0324ebeSBarry Smith     int level=3;
1113*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_StructPFMGSetPrintLevel(ex->hsolver,level));
111468326731SBarry Smith   }
11159e5bc791SBarry 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);
1116*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGSetMaxIter(ex->hsolver,ex->its));
11179e5bc791SBarry 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);
1118*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGSetNumPreRelax(ex->hsolver,ex->num_pre_relax));
11199e5bc791SBarry 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);
1120*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGSetNumPostRelax(ex->hsolver,ex->num_post_relax));
11219e5bc791SBarry Smith 
11223b46a515SGlenn Hammond   ierr = PetscOptionsInt("-pc_pfmg_max_levels","Max Levels for MG hierarchy","HYPRE_StructPFMGSetMaxLevels",ex->max_levels,&ex->max_levels,PETSC_NULL);CHKERRQ(ierr);
1123*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGSetMaxLevels(ex->hsolver,ex->max_levels));
11243b46a515SGlenn Hammond 
11259e5bc791SBarry Smith   ierr = PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr);
1126*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGSetTol(ex->hsolver,ex->tol));
11279e5bc791SBarry 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);
1128*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGSetRelaxType(ex->hsolver, ex->relax_type));
11299e5bc791SBarry Smith   ierr = PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,2,PFMGRAPType[ex->rap_type],&ex->rap_type,PETSC_NULL);CHKERRQ(ierr);
1130*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGSetRAPType(ex->hsolver, ex->rap_type));
1131ebc551c0SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1132ebc551c0SBarry Smith   PetscFunctionReturn(0);
1133ebc551c0SBarry Smith }
1134ebc551c0SBarry Smith 
1135f91d8e95SBarry Smith #undef __FUNCT__
1136f91d8e95SBarry Smith #define __FUNCT__ "PCApply_PFMG"
1137f91d8e95SBarry Smith PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y)
1138f91d8e95SBarry Smith {
1139f91d8e95SBarry Smith   PetscErrorCode  ierr;
1140f91d8e95SBarry Smith   PC_PFMG         *ex = (PC_PFMG*) pc->data;
1141f91d8e95SBarry Smith   PetscScalar     *xx,*yy;
1142f91d8e95SBarry Smith   int             ilower[3],iupper[3];
114368326731SBarry Smith   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
1144f91d8e95SBarry Smith 
1145f91d8e95SBarry Smith   PetscFunctionBegin;
114668326731SBarry Smith   ierr = DAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1147f91d8e95SBarry Smith   iupper[0] += ilower[0] - 1;
1148f91d8e95SBarry Smith   iupper[1] += ilower[1] - 1;
1149f91d8e95SBarry Smith   iupper[2] += ilower[2] - 1;
1150f91d8e95SBarry Smith 
1151f91d8e95SBarry Smith   /* copy x values over to hypre */
1152*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructVectorSetConstantValues(mx->hb,0.0));
1153f91d8e95SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1154*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructVectorSetBoxValues(mx->hb,ilower,iupper,xx));
1155f91d8e95SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1156*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructVectorAssemble(mx->hb));
1157*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGSolve(ex->hsolver,mx->hmat,mx->hb,mx->hx));
1158f91d8e95SBarry Smith 
1159f91d8e95SBarry Smith   /* copy solution values back to PETSc */
1160f91d8e95SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1161*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructVectorGetBoxValues(mx->hx,ilower,iupper,yy));
1162f91d8e95SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1163f91d8e95SBarry Smith   PetscFunctionReturn(0);
1164f91d8e95SBarry Smith }
1165f91d8e95SBarry Smith 
11669e5bc791SBarry Smith #undef __FUNCT__
11679e5bc791SBarry Smith #define __FUNCT__ "PCApplyRichardson_PFMG"
11687319c654SBarry Smith static PetscErrorCode PCApplyRichardson_PFMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
11699e5bc791SBarry Smith {
11709e5bc791SBarry Smith   PC_PFMG        *jac = (PC_PFMG*)pc->data;
11719e5bc791SBarry Smith   PetscErrorCode ierr;
11729e5bc791SBarry Smith   int            oits;
11739e5bc791SBarry Smith 
11749e5bc791SBarry Smith   PetscFunctionBegin;
1175*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGSetMaxIter(jac->hsolver,its*jac->its));
1176*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGSetTol(jac->hsolver,rtol));
11779e5bc791SBarry Smith 
11789e5bc791SBarry Smith   ierr = PCApply_PFMG(pc,b,y);CHKERRQ(ierr);
1179*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGGetNumIterations(jac->hsolver,&oits));
11809e5bc791SBarry Smith   *outits = oits;
11819e5bc791SBarry Smith   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
11829e5bc791SBarry Smith   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
1183*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGSetTol(jac->hsolver,jac->tol));
1184*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGSetMaxIter(jac->hsolver,jac->its));
11859e5bc791SBarry Smith   PetscFunctionReturn(0);
11869e5bc791SBarry Smith }
11879e5bc791SBarry Smith 
11889e5bc791SBarry Smith 
11893a32d3dbSGlenn Hammond #undef __FUNCT__
11903a32d3dbSGlenn Hammond #define __FUNCT__ "PCSetUp_PFMG"
11913a32d3dbSGlenn Hammond PetscErrorCode PCSetUp_PFMG(PC pc)
11923a32d3dbSGlenn Hammond {
11933a32d3dbSGlenn Hammond   PetscErrorCode  ierr;
11943a32d3dbSGlenn Hammond   PC_PFMG         *ex = (PC_PFMG*) pc->data;
11953a32d3dbSGlenn Hammond   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(pc->pmat->data);
11963a32d3dbSGlenn Hammond   PetscTruth      flg;
11973a32d3dbSGlenn Hammond 
11983a32d3dbSGlenn Hammond   PetscFunctionBegin;
11993a32d3dbSGlenn Hammond   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATHYPRESTRUCT,&flg);CHKERRQ(ierr);
1200e7e72b3dSBarry Smith   if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESTRUCT with this preconditioner");
12013a32d3dbSGlenn Hammond 
12023a32d3dbSGlenn Hammond   /* create the hypre solver object and set its information */
12033a32d3dbSGlenn Hammond   if (ex->hsolver) {
1204*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_StructPFMGDestroy(ex->hsolver));
12053a32d3dbSGlenn Hammond   }
1206*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGCreate(ex->hcomm,&ex->hsolver));
12073a32d3dbSGlenn Hammond   ierr = PCSetFromOptions_PFMG(pc);CHKERRQ(ierr);
1208*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGSetup(ex->hsolver,mx->hmat,mx->hb,mx->hx));
1209*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGSetZeroGuess(ex->hsolver));
12103a32d3dbSGlenn Hammond   PetscFunctionReturn(0);
12113a32d3dbSGlenn Hammond }
12123a32d3dbSGlenn Hammond 
1213ebc551c0SBarry Smith 
1214ebc551c0SBarry Smith /*MC
1215ebc551c0SBarry Smith      PCPFMG - the hypre PFMG multigrid solver
1216ebc551c0SBarry Smith 
1217ebc551c0SBarry Smith    Level: advanced
1218ebc551c0SBarry Smith 
12199e5bc791SBarry Smith    Options Database:
12209e5bc791SBarry Smith + -pc_pfmg_its <its> number of iterations of PFMG to use as preconditioner
12219e5bc791SBarry Smith . -pc_pfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
12229e5bc791SBarry Smith . -pc_pfmg_num_post_relax <steps> number of smoothing steps after coarse grid
12239e5bc791SBarry Smith . -pc_pfmg_tol <tol> tolerance of PFMG
12249e5bc791SBarry 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
12259e5bc791SBarry Smith - -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin
1226f91d8e95SBarry Smith 
12279e5bc791SBarry Smith    Notes:  This is for CELL-centered descretizations
12289e5bc791SBarry Smith 
12298e395302SJed Brown            This must be used with the MATHYPRESTRUCT matrix type.
12309e5bc791SBarry Smith            This is less general than in hypre, it supports only one block per process defined by a PETSc DA.
12319e5bc791SBarry Smith 
12329e5bc791SBarry Smith .seealso:  PCMG, MATHYPRESTRUCT
1233ebc551c0SBarry Smith M*/
1234ebc551c0SBarry Smith 
1235ebc551c0SBarry Smith EXTERN_C_BEGIN
1236ebc551c0SBarry Smith #undef __FUNCT__
1237ebc551c0SBarry Smith #define __FUNCT__ "PCCreate_PFMG"
1238ebc551c0SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_PFMG(PC pc)
1239ebc551c0SBarry Smith {
1240ebc551c0SBarry Smith   PetscErrorCode ierr;
1241ebc551c0SBarry Smith   PC_PFMG        *ex;
1242ebc551c0SBarry Smith 
1243ebc551c0SBarry Smith   PetscFunctionBegin;
1244ebc551c0SBarry Smith   ierr = PetscNew(PC_PFMG,&ex);CHKERRQ(ierr);\
124568326731SBarry Smith   pc->data = ex;
1246ebc551c0SBarry Smith 
12479e5bc791SBarry Smith   ex->its            = 1;
12489e5bc791SBarry Smith   ex->tol            = 1.e-8;
12499e5bc791SBarry Smith   ex->relax_type     = 1;
12509e5bc791SBarry Smith   ex->rap_type       = 0;
12519e5bc791SBarry Smith   ex->num_pre_relax  = 1;
12529e5bc791SBarry Smith   ex->num_post_relax = 1;
12533b46a515SGlenn Hammond   ex->max_levels     = 0;
12549e5bc791SBarry Smith 
1255ebc551c0SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_PFMG;
1256ebc551c0SBarry Smith   pc->ops->view            = PCView_PFMG;
1257ebc551c0SBarry Smith   pc->ops->destroy         = PCDestroy_PFMG;
1258f91d8e95SBarry Smith   pc->ops->apply           = PCApply_PFMG;
12599e5bc791SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_PFMG;
126068326731SBarry Smith   pc->ops->setup           = PCSetUp_PFMG;
1261f91d8e95SBarry Smith   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr);
1262*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_StructPFMGCreate(ex->hcomm,&ex->hsolver));
1263ebc551c0SBarry Smith   PetscFunctionReturn(0);
1264ebc551c0SBarry Smith }
1265ebc551c0SBarry Smith EXTERN_C_END
1266d851a50bSGlenn Hammond 
1267d851a50bSGlenn Hammond /* we know we are working with a HYPRE_SStructMatrix */
1268d851a50bSGlenn Hammond typedef struct {
1269d851a50bSGlenn Hammond   MPI_Comm            hcomm;       /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */
1270d851a50bSGlenn Hammond   HYPRE_SStructSolver  ss_solver;
1271d851a50bSGlenn Hammond 
1272d851a50bSGlenn Hammond   /* keep copy of SYSPFMG options used so may view them */
1273d851a50bSGlenn Hammond   int                 its;
1274d851a50bSGlenn Hammond   double              tol;
1275d851a50bSGlenn Hammond   int                 relax_type;
1276d851a50bSGlenn Hammond   int                 num_pre_relax,num_post_relax;
1277d851a50bSGlenn Hammond } PC_SysPFMG;
1278d851a50bSGlenn Hammond 
1279d851a50bSGlenn Hammond #undef __FUNCT__
1280d851a50bSGlenn Hammond #define __FUNCT__ "PCDestroy_SysPFMG"
1281d851a50bSGlenn Hammond PetscErrorCode PCDestroy_SysPFMG(PC pc)
1282d851a50bSGlenn Hammond {
1283d851a50bSGlenn Hammond   PetscErrorCode ierr;
1284d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1285d851a50bSGlenn Hammond 
1286d851a50bSGlenn Hammond   PetscFunctionBegin;
1287*a65764d7SBarry Smith   if (ex->ss_solver) {PetscStackCallHypre("",HYPRE_SStructSysPFMGDestroy(ex->ss_solver));}
1288d851a50bSGlenn Hammond   ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr);
1289d851a50bSGlenn Hammond   ierr = PetscFree(ex);CHKERRQ(ierr);
1290d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1291d851a50bSGlenn Hammond }
1292d851a50bSGlenn Hammond 
1293d851a50bSGlenn Hammond static const char *SysPFMGRelaxType[]   = {"Weighted-Jacobi","Red/Black-Gauss-Seidel"};
1294d851a50bSGlenn Hammond 
1295d851a50bSGlenn Hammond #undef __FUNCT__
1296d851a50bSGlenn Hammond #define __FUNCT__ "PCView_SysPFMG"
1297d851a50bSGlenn Hammond PetscErrorCode PCView_SysPFMG(PC pc,PetscViewer viewer)
1298d851a50bSGlenn Hammond {
1299d851a50bSGlenn Hammond   PetscErrorCode ierr;
1300d851a50bSGlenn Hammond   PetscTruth     iascii;
1301d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1302d851a50bSGlenn Hammond 
1303d851a50bSGlenn Hammond   PetscFunctionBegin;
13042692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1305d851a50bSGlenn Hammond   if (iascii) {
1306d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG preconditioning\n");CHKERRQ(ierr);
1307d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: max iterations %d\n",ex->its);CHKERRQ(ierr);
1308d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr);
1309d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr);
1310d851a50bSGlenn Hammond     ierr = PetscViewerASCIIPrintf(viewer,"  HYPRE SysPFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr);
1311d851a50bSGlenn Hammond   }
1312d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1313d851a50bSGlenn Hammond }
1314d851a50bSGlenn Hammond 
1315d851a50bSGlenn Hammond 
1316d851a50bSGlenn Hammond #undef __FUNCT__
1317d851a50bSGlenn Hammond #define __FUNCT__ "PCSetFromOptions_SysPFMG"
1318d851a50bSGlenn Hammond PetscErrorCode PCSetFromOptions_SysPFMG(PC pc)
1319d851a50bSGlenn Hammond {
1320d851a50bSGlenn Hammond   PetscErrorCode ierr;
1321d851a50bSGlenn Hammond   PC_SysPFMG    *ex = (PC_SysPFMG*) pc->data;
1322d851a50bSGlenn Hammond   PetscTruth     flg = PETSC_FALSE;
1323d851a50bSGlenn Hammond 
1324d851a50bSGlenn Hammond   PetscFunctionBegin;
1325d851a50bSGlenn Hammond   ierr = PetscOptionsHead("SysPFMG options");CHKERRQ(ierr);
1326d851a50bSGlenn Hammond   ierr = PetscOptionsTruth("-pc_syspfmg_print_statistics","Print statistics","HYPRE_SStructSysPFMGSetPrintLevel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1327d851a50bSGlenn Hammond   if (flg) {
1328d851a50bSGlenn Hammond     int level=3;
1329*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_SStructSysPFMGSetPrintLevel(ex->ss_solver,level));
1330d851a50bSGlenn Hammond   }
1331d851a50bSGlenn 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);
1332*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_SStructSysPFMGSetMaxIter(ex->ss_solver,ex->its));
1333d851a50bSGlenn 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);
1334*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_SStructSysPFMGSetNumPreRelax(ex->ss_solver,ex->num_pre_relax));
1335d851a50bSGlenn 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);
1336*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_SStructSysPFMGSetNumPostRelax(ex->ss_solver,ex->num_post_relax));
1337d851a50bSGlenn Hammond 
1338d851a50bSGlenn Hammond   ierr = PetscOptionsReal("-pc_syspfmg_tol","Tolerance of SysPFMG","HYPRE_SStructSysPFMGSetTol",ex->tol,&ex->tol,PETSC_NULL);CHKERRQ(ierr);
1339*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_SStructSysPFMGSetTol(ex->ss_solver,ex->tol));
1340d851a50bSGlenn 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);
1341*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_SStructSysPFMGSetRelaxType(ex->ss_solver, ex->relax_type));
1342d851a50bSGlenn Hammond   ierr = PetscOptionsTail();CHKERRQ(ierr);
1343d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1344d851a50bSGlenn Hammond }
1345d851a50bSGlenn Hammond 
1346d851a50bSGlenn Hammond #undef __FUNCT__
1347d851a50bSGlenn Hammond #define __FUNCT__ "PCApply_SysPFMG"
1348d851a50bSGlenn Hammond PetscErrorCode PCApply_SysPFMG(PC pc,Vec x,Vec y)
1349d851a50bSGlenn Hammond {
1350d851a50bSGlenn Hammond   PetscErrorCode    ierr;
1351d851a50bSGlenn Hammond   PC_SysPFMG       *ex = (PC_SysPFMG*) pc->data;
1352d851a50bSGlenn Hammond   PetscScalar      *xx,*yy;
1353d851a50bSGlenn Hammond   int               ilower[3],iupper[3];
1354d851a50bSGlenn Hammond   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct *)(pc->pmat->data);
1355d851a50bSGlenn Hammond   int               ordering= mx->dofs_order;
1356d851a50bSGlenn Hammond   int               nvars= mx->nvars;
1357d851a50bSGlenn Hammond   int               part= 0;
1358d851a50bSGlenn Hammond   int               size;
1359d851a50bSGlenn Hammond   int               i;
1360d851a50bSGlenn Hammond 
1361d851a50bSGlenn Hammond   PetscFunctionBegin;
1362d851a50bSGlenn Hammond   ierr = DAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1363d851a50bSGlenn Hammond   iupper[0] += ilower[0] - 1;
1364d851a50bSGlenn Hammond   iupper[1] += ilower[1] - 1;
1365d851a50bSGlenn Hammond   iupper[2] += ilower[2] - 1;
1366d851a50bSGlenn Hammond 
1367d851a50bSGlenn Hammond   size= 1;
1368d851a50bSGlenn Hammond   for (i= 0; i< 3; i++) {
1369d851a50bSGlenn Hammond      size*= (iupper[i]-ilower[i]+1);
1370d851a50bSGlenn Hammond   }
1371d851a50bSGlenn Hammond   /* copy x values over to hypre for variable ordering */
1372d851a50bSGlenn Hammond   if (ordering) {
1373*a65764d7SBarry Smith      PetscStackCallHypre("",HYPRE_SStructVectorSetConstantValues(mx->ss_b,0.0));
1374d851a50bSGlenn Hammond      ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1375d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
1376*a65764d7SBarry Smith         PetscStackCallHypre("",HYPRE_SStructVectorSetBoxValues(mx->ss_b,part,ilower,iupper,i,xx+(size*i)));
1377d851a50bSGlenn Hammond      }
1378d851a50bSGlenn Hammond      ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1379*a65764d7SBarry Smith      PetscStackCallHypre("",HYPRE_SStructVectorAssemble(mx->ss_b));
1380*a65764d7SBarry Smith      PetscStackCallHypre("",HYPRE_SStructMatrixMatvec(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
1381*a65764d7SBarry Smith      PetscStackCallHypre("",HYPRE_SStructSysPFMGSolve(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1382d851a50bSGlenn Hammond 
1383d851a50bSGlenn Hammond      /* copy solution values back to PETSc */
1384d851a50bSGlenn Hammond      ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1385d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
1386*a65764d7SBarry Smith         PetscStackCallHypre("",HYPRE_SStructVectorGetBoxValues(mx->ss_x,part,ilower,iupper,i,yy+(size*i)));
1387d851a50bSGlenn Hammond      }
1388d851a50bSGlenn Hammond      ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1389*a65764d7SBarry Smith   } else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
1390d851a50bSGlenn Hammond      PetscScalar     *z;
1391d851a50bSGlenn Hammond      int              j, k;
1392d851a50bSGlenn Hammond 
1393d851a50bSGlenn Hammond      ierr = PetscMalloc(nvars*size*sizeof(PetscScalar),&z);CHKERRQ(ierr);
1394*a65764d7SBarry Smith      PetscStackCallHypre("",HYPRE_SStructVectorSetConstantValues(mx->ss_b,0.0));
1395d851a50bSGlenn Hammond      ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1396d851a50bSGlenn Hammond 
1397d851a50bSGlenn Hammond      /* transform nodal to hypre's variable ordering for sys_pfmg */
1398d851a50bSGlenn Hammond      for (i= 0; i< size; i++) {
1399d851a50bSGlenn Hammond         k= i*nvars;
1400d851a50bSGlenn Hammond         for (j= 0; j< nvars; j++) {
1401d851a50bSGlenn Hammond            z[j*size+i]= xx[k+j];
1402d851a50bSGlenn Hammond         }
1403d851a50bSGlenn Hammond      }
1404d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
1405*a65764d7SBarry Smith         PetscStackCallHypre("",HYPRE_SStructVectorSetBoxValues(mx->ss_b,part,ilower,iupper,i,z+(size*i)));
1406d851a50bSGlenn Hammond      }
1407d851a50bSGlenn Hammond      ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1408*a65764d7SBarry Smith      PetscStackCallHypre("",HYPRE_SStructVectorAssemble(mx->ss_b));
1409*a65764d7SBarry Smith      PetscStackCallHypre("",HYPRE_SStructSysPFMGSolve(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1410d851a50bSGlenn Hammond 
1411d851a50bSGlenn Hammond      /* copy solution values back to PETSc */
1412d851a50bSGlenn Hammond      ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1413d851a50bSGlenn Hammond      for (i= 0; i< nvars; i++) {
1414*a65764d7SBarry Smith         PetscStackCallHypre("",HYPRE_SStructVectorGetBoxValues(mx->ss_x,part,ilower,iupper,i,z+(size*i)));
1415d851a50bSGlenn Hammond      }
1416d851a50bSGlenn Hammond      /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
1417d851a50bSGlenn Hammond      for (i= 0; i< size; i++) {
1418d851a50bSGlenn Hammond         k= i*nvars;
1419d851a50bSGlenn Hammond         for (j= 0; j< nvars; j++) {
1420d851a50bSGlenn Hammond            yy[k+j]= z[j*size+i];
1421d851a50bSGlenn Hammond         }
1422d851a50bSGlenn Hammond      }
1423d851a50bSGlenn Hammond      ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1424d851a50bSGlenn Hammond      ierr = PetscFree(z);CHKERRQ(ierr);
1425d851a50bSGlenn Hammond   }
1426d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1427d851a50bSGlenn Hammond }
1428d851a50bSGlenn Hammond 
1429d851a50bSGlenn Hammond #undef __FUNCT__
1430d851a50bSGlenn Hammond #define __FUNCT__ "PCApplyRichardson_SysPFMG"
1431d851a50bSGlenn Hammond static PetscErrorCode PCApplyRichardson_SysPFMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason)
1432d851a50bSGlenn Hammond {
1433d851a50bSGlenn Hammond   PC_SysPFMG    *jac = (PC_SysPFMG*)pc->data;
1434d851a50bSGlenn Hammond   PetscErrorCode ierr;
1435d851a50bSGlenn Hammond   int            oits;
1436d851a50bSGlenn Hammond 
1437d851a50bSGlenn Hammond   PetscFunctionBegin;
1438*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_SStructSysPFMGSetMaxIter(jac->ss_solver,its*jac->its));
1439*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_SStructSysPFMGSetTol(jac->ss_solver,rtol));
1440d851a50bSGlenn Hammond   ierr = PCApply_SysPFMG(pc,b,y);CHKERRQ(ierr);
1441*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_SStructSysPFMGGetNumIterations(jac->ss_solver,&oits));
1442d851a50bSGlenn Hammond   *outits = oits;
1443d851a50bSGlenn Hammond   if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS;
1444d851a50bSGlenn Hammond   else             *reason = PCRICHARDSON_CONVERGED_RTOL;
1445*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_SStructSysPFMGSetTol(jac->ss_solver,jac->tol));
1446*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_SStructSysPFMGSetMaxIter(jac->ss_solver,jac->its));
1447d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1448d851a50bSGlenn Hammond }
1449d851a50bSGlenn Hammond 
1450d851a50bSGlenn Hammond 
1451d851a50bSGlenn Hammond #undef __FUNCT__
1452d851a50bSGlenn Hammond #define __FUNCT__ "PCSetUp_SysPFMG"
1453d851a50bSGlenn Hammond PetscErrorCode PCSetUp_SysPFMG(PC pc)
1454d851a50bSGlenn Hammond {
1455d851a50bSGlenn Hammond   PetscErrorCode    ierr;
1456d851a50bSGlenn Hammond   PC_SysPFMG        *ex = (PC_SysPFMG*) pc->data;
1457d851a50bSGlenn Hammond   Mat_HYPRESStruct  *mx = (Mat_HYPRESStruct *)(pc->pmat->data);
1458d851a50bSGlenn Hammond   PetscTruth        flg;
1459d851a50bSGlenn Hammond 
1460d851a50bSGlenn Hammond   PetscFunctionBegin;
1461d851a50bSGlenn Hammond   ierr = PetscTypeCompare((PetscObject)pc->pmat,MATHYPRESSTRUCT,&flg);CHKERRQ(ierr);
1462e7e72b3dSBarry Smith   if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESSTRUCT with this preconditioner");
1463d851a50bSGlenn Hammond 
1464d851a50bSGlenn Hammond   /* create the hypre sstruct solver object and set its information */
1465d851a50bSGlenn Hammond   if (ex->ss_solver) {
1466*a65764d7SBarry Smith     PetscStackCallHypre("",HYPRE_SStructSysPFMGDestroy(ex->ss_solver));
1467d851a50bSGlenn Hammond   }
1468*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_SStructSysPFMGCreate(ex->hcomm,&ex->ss_solver));
1469d851a50bSGlenn Hammond   ierr = PCSetFromOptions_SysPFMG(pc);CHKERRQ(ierr);
1470*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_SStructSysPFMGSetZeroGuess(ex->ss_solver));
1471*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_SStructSysPFMGSetup(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x));
1472d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1473d851a50bSGlenn Hammond }
1474d851a50bSGlenn Hammond 
1475d851a50bSGlenn Hammond 
1476d851a50bSGlenn Hammond /*MC
1477d851a50bSGlenn Hammond      PCSysPFMG - the hypre SysPFMG multigrid solver
1478d851a50bSGlenn Hammond 
1479d851a50bSGlenn Hammond    Level: advanced
1480d851a50bSGlenn Hammond 
1481d851a50bSGlenn Hammond    Options Database:
1482d851a50bSGlenn Hammond + -pc_syspfmg_its <its> number of iterations of SysPFMG to use as preconditioner
1483d851a50bSGlenn Hammond . -pc_syspfmg_num_pre_relax <steps> number of smoothing steps before coarse grid
1484d851a50bSGlenn Hammond . -pc_syspfmg_num_post_relax <steps> number of smoothing steps after coarse grid
1485d851a50bSGlenn Hammond . -pc_syspfmg_tol <tol> tolerance of SysPFMG
1486d851a50bSGlenn Hammond . -pc_syspfmg_relax_type -relaxation type for the up and down cycles, one of Weighted-Jacobi,Red/Black-Gauss-Seidel
1487d851a50bSGlenn Hammond 
1488d851a50bSGlenn Hammond    Notes:  This is for CELL-centered descretizations
1489d851a50bSGlenn Hammond 
1490f6680f47SSatish Balay            This must be used with the MATHYPRESSTRUCT matrix type.
1491d851a50bSGlenn Hammond            This is less general than in hypre, it supports only one part, and one block per process defined by a PETSc DA.
1492d851a50bSGlenn Hammond            Also, only cell-centered variables.
1493d851a50bSGlenn Hammond 
1494d851a50bSGlenn Hammond .seealso:  PCMG, MATHYPRESSTRUCT
1495d851a50bSGlenn Hammond M*/
1496d851a50bSGlenn Hammond 
1497d851a50bSGlenn Hammond EXTERN_C_BEGIN
1498d851a50bSGlenn Hammond #undef __FUNCT__
1499d851a50bSGlenn Hammond #define __FUNCT__ "PCCreate_SysPFMG"
1500d851a50bSGlenn Hammond PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_SysPFMG(PC pc)
1501d851a50bSGlenn Hammond {
1502d851a50bSGlenn Hammond   PetscErrorCode     ierr;
1503d851a50bSGlenn Hammond   PC_SysPFMG        *ex;
1504d851a50bSGlenn Hammond 
1505d851a50bSGlenn Hammond   PetscFunctionBegin;
1506d851a50bSGlenn Hammond   ierr = PetscNew(PC_SysPFMG,&ex);CHKERRQ(ierr);\
1507d851a50bSGlenn Hammond   pc->data = ex;
1508d851a50bSGlenn Hammond 
1509d851a50bSGlenn Hammond   ex->its            = 1;
1510d851a50bSGlenn Hammond   ex->tol            = 1.e-8;
1511d851a50bSGlenn Hammond   ex->relax_type     = 1;
1512d851a50bSGlenn Hammond   ex->num_pre_relax  = 1;
1513d851a50bSGlenn Hammond   ex->num_post_relax = 1;
1514d851a50bSGlenn Hammond 
1515d851a50bSGlenn Hammond   pc->ops->setfromoptions  = PCSetFromOptions_SysPFMG;
1516d851a50bSGlenn Hammond   pc->ops->view            = PCView_SysPFMG;
1517d851a50bSGlenn Hammond   pc->ops->destroy         = PCDestroy_SysPFMG;
1518d851a50bSGlenn Hammond   pc->ops->apply           = PCApply_SysPFMG;
1519d851a50bSGlenn Hammond   pc->ops->applyrichardson = PCApplyRichardson_SysPFMG;
1520d851a50bSGlenn Hammond   pc->ops->setup           = PCSetUp_SysPFMG;
1521d851a50bSGlenn Hammond   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ex->hcomm));CHKERRQ(ierr);
1522*a65764d7SBarry Smith   PetscStackCallHypre("",HYPRE_SStructSysPFMGCreate(ex->hcomm,&ex->ss_solver));
1523d851a50bSGlenn Hammond   PetscFunctionReturn(0);
1524d851a50bSGlenn Hammond }
1525d851a50bSGlenn Hammond EXTERN_C_END
1526