1*16d9e3a6SLisandro Dalcin #define PETSCKSP_DLL 2*16d9e3a6SLisandro Dalcin 3*16d9e3a6SLisandro Dalcin /* 4*16d9e3a6SLisandro Dalcin Provides an interface to the LLNL package hypre 5*16d9e3a6SLisandro Dalcin */ 6*16d9e3a6SLisandro Dalcin #include "private/pcimpl.h" /*I "petscpc.h" I*/ 7*16d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 8*16d9e3a6SLisandro Dalcin #include "HYPRE.h" 9*16d9e3a6SLisandro Dalcin #include "IJ_mv.h" 10*16d9e3a6SLisandro Dalcin #include "parcsr_ls.h" 11*16d9e3a6SLisandro Dalcin EXTERN_C_END 12*16d9e3a6SLisandro Dalcin 13*16d9e3a6SLisandro Dalcin EXTERN PetscErrorCode MatHYPRE_IJMatrixCreate(Mat,HYPRE_IJMatrix*); 14*16d9e3a6SLisandro Dalcin EXTERN PetscErrorCode MatHYPRE_IJMatrixCopy(Mat,HYPRE_IJMatrix); 15*16d9e3a6SLisandro Dalcin EXTERN PetscErrorCode VecHYPRE_IJVectorCreate(Vec,HYPRE_IJVector*); 16*16d9e3a6SLisandro Dalcin 17*16d9e3a6SLisandro Dalcin /* 18*16d9e3a6SLisandro Dalcin Private context (data structure) for the preconditioner. 19*16d9e3a6SLisandro Dalcin */ 20*16d9e3a6SLisandro Dalcin typedef struct { 21*16d9e3a6SLisandro Dalcin HYPRE_Solver hsolver; 22*16d9e3a6SLisandro Dalcin HYPRE_IJMatrix ij; 23*16d9e3a6SLisandro Dalcin HYPRE_IJVector b,x; 24*16d9e3a6SLisandro Dalcin 25*16d9e3a6SLisandro Dalcin PetscErrorCode (*destroy)(HYPRE_Solver); 26*16d9e3a6SLisandro Dalcin PetscErrorCode (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector); 27*16d9e3a6SLisandro Dalcin PetscErrorCode (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector); 28*16d9e3a6SLisandro Dalcin 29*16d9e3a6SLisandro Dalcin MPI_Comm comm_hypre; 30*16d9e3a6SLisandro Dalcin char *hypre_type; 31*16d9e3a6SLisandro Dalcin 32*16d9e3a6SLisandro Dalcin /* options for Pilut and BoomerAMG*/ 33*16d9e3a6SLisandro Dalcin int maxiter; 34*16d9e3a6SLisandro Dalcin double tol; 35*16d9e3a6SLisandro Dalcin 36*16d9e3a6SLisandro Dalcin /* options for Pilut */ 37*16d9e3a6SLisandro Dalcin int factorrowsize; 38*16d9e3a6SLisandro Dalcin 39*16d9e3a6SLisandro Dalcin /* options for ParaSails */ 40*16d9e3a6SLisandro Dalcin int nlevels; 41*16d9e3a6SLisandro Dalcin double threshhold; 42*16d9e3a6SLisandro Dalcin double filter; 43*16d9e3a6SLisandro Dalcin int sym; 44*16d9e3a6SLisandro Dalcin double loadbal; 45*16d9e3a6SLisandro Dalcin int logging; 46*16d9e3a6SLisandro Dalcin int ruse; 47*16d9e3a6SLisandro Dalcin int symt; 48*16d9e3a6SLisandro Dalcin 49*16d9e3a6SLisandro Dalcin /* options for Euclid */ 50*16d9e3a6SLisandro Dalcin PetscTruth bjilu; 51*16d9e3a6SLisandro Dalcin int levels; 52*16d9e3a6SLisandro Dalcin 53*16d9e3a6SLisandro Dalcin /* options for Euclid and BoomerAMG */ 54*16d9e3a6SLisandro Dalcin PetscTruth printstatistics; 55*16d9e3a6SLisandro Dalcin 56*16d9e3a6SLisandro Dalcin /* options for BoomerAMG */ 57*16d9e3a6SLisandro Dalcin int cycletype; 58*16d9e3a6SLisandro Dalcin int maxlevels; 59*16d9e3a6SLisandro Dalcin double strongthreshold; 60*16d9e3a6SLisandro Dalcin double maxrowsum; 61*16d9e3a6SLisandro Dalcin int gridsweeps[4]; 62*16d9e3a6SLisandro Dalcin int coarsentype; 63*16d9e3a6SLisandro Dalcin int measuretype; 64*16d9e3a6SLisandro Dalcin int relaxtype[4]; 65*16d9e3a6SLisandro Dalcin double relaxweight; 66*16d9e3a6SLisandro Dalcin double outerrelaxweight; 67*16d9e3a6SLisandro Dalcin int relaxorder; 68*16d9e3a6SLisandro Dalcin double truncfactor; 69*16d9e3a6SLisandro Dalcin PetscTruth applyrichardson; 70*16d9e3a6SLisandro Dalcin 71*16d9e3a6SLisandro Dalcin } PC_HYPRE; 72*16d9e3a6SLisandro Dalcin 73*16d9e3a6SLisandro Dalcin 74*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 75*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetUp_HYPRE" 76*16d9e3a6SLisandro Dalcin static PetscErrorCode PCSetUp_HYPRE(PC pc) 77*16d9e3a6SLisandro Dalcin { 78*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 79*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 80*16d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 81*16d9e3a6SLisandro Dalcin HYPRE_ParVector bv,xv; 82*16d9e3a6SLisandro Dalcin PetscInt bs; 83*16d9e3a6SLisandro Dalcin int hierr; 84*16d9e3a6SLisandro Dalcin 85*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 86*16d9e3a6SLisandro Dalcin if (!jac->hypre_type) { 87*16d9e3a6SLisandro Dalcin ierr = PCHYPRESetType(pc,"pilut");CHKERRQ(ierr); 88*16d9e3a6SLisandro Dalcin } 89*16d9e3a6SLisandro Dalcin #if 1 90*16d9e3a6SLisandro Dalcin if (!pc->setupcalled) { 91*16d9e3a6SLisandro Dalcin /* create the matrix and vectors the first time through */ 92*16d9e3a6SLisandro Dalcin Vec x,b; 93*16d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr); 94*16d9e3a6SLisandro Dalcin ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr); 95*16d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr); 96*16d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr); 97*16d9e3a6SLisandro Dalcin ierr = VecDestroy(x);CHKERRQ(ierr); 98*16d9e3a6SLisandro Dalcin ierr = VecDestroy(b);CHKERRQ(ierr); 99*16d9e3a6SLisandro Dalcin } else if (pc->flag != SAME_NONZERO_PATTERN) { 100*16d9e3a6SLisandro Dalcin /* rebuild the matrix from scratch */ 101*16d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); 102*16d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr); 103*16d9e3a6SLisandro Dalcin } 104*16d9e3a6SLisandro Dalcin #else 105*16d9e3a6SLisandro Dalcin if (!jac->ij) { /* create the matrix the first time through */ 106*16d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr); 107*16d9e3a6SLisandro Dalcin } 108*16d9e3a6SLisandro Dalcin if (!jac->b) { /* create the vectors the first time through */ 109*16d9e3a6SLisandro Dalcin Vec x,b; 110*16d9e3a6SLisandro Dalcin ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr); 111*16d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr); 112*16d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr); 113*16d9e3a6SLisandro Dalcin ierr = VecDestroy(x);CHKERRQ(ierr); 114*16d9e3a6SLisandro Dalcin ierr = VecDestroy(b);CHKERRQ(ierr); 115*16d9e3a6SLisandro Dalcin } 116*16d9e3a6SLisandro Dalcin #endif 117*16d9e3a6SLisandro Dalcin /* special case for BoomerAMG */ 118*16d9e3a6SLisandro Dalcin if (jac->setup == HYPRE_BoomerAMGSetup) { 119*16d9e3a6SLisandro Dalcin ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr); 120*16d9e3a6SLisandro Dalcin if (bs > 1) { 121*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetNumFunctions(jac->hsolver,bs);CHKERRQ(ierr); 122*16d9e3a6SLisandro Dalcin } 123*16d9e3a6SLisandro Dalcin }; 124*16d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr); 125*16d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 126*16d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&bv);CHKERRQ(ierr); 127*16d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&xv);CHKERRQ(ierr); 128*16d9e3a6SLisandro Dalcin hierr = (*jac->setup)(jac->hsolver,hmat,bv,xv); 129*16d9e3a6SLisandro Dalcin if (hierr) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE setup, error code %d",hierr); 130*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 131*16d9e3a6SLisandro Dalcin } 132*16d9e3a6SLisandro Dalcin 133*16d9e3a6SLisandro Dalcin /* 134*16d9e3a6SLisandro Dalcin Replaces the address where the HYPRE vector points to its data with the address of 135*16d9e3a6SLisandro Dalcin PETSc's data. Saves the old address so it can be reset when we are finished with it. 136*16d9e3a6SLisandro Dalcin Allows use to get the data into a HYPRE vector without the cost of memcopies 137*16d9e3a6SLisandro Dalcin */ 138*16d9e3a6SLisandro Dalcin #define HYPREReplacePointer(b,newvalue,savedvalue) {\ 139*16d9e3a6SLisandro Dalcin hypre_ParVector *par_vector = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));\ 140*16d9e3a6SLisandro Dalcin hypre_Vector *local_vector = hypre_ParVectorLocalVector(par_vector);\ 141*16d9e3a6SLisandro Dalcin savedvalue = local_vector->data;\ 142*16d9e3a6SLisandro Dalcin local_vector->data = newvalue;} 143*16d9e3a6SLisandro Dalcin 144*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 145*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApply_HYPRE" 146*16d9e3a6SLisandro Dalcin static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x) 147*16d9e3a6SLisandro Dalcin { 148*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 149*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 150*16d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 151*16d9e3a6SLisandro Dalcin PetscScalar *bv,*xv; 152*16d9e3a6SLisandro Dalcin HYPRE_ParVector jbv,jxv; 153*16d9e3a6SLisandro Dalcin PetscScalar *sbv,*sxv; 154*16d9e3a6SLisandro Dalcin int hierr; 155*16d9e3a6SLisandro Dalcin 156*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 157*16d9e3a6SLisandro Dalcin if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);} 158*16d9e3a6SLisandro Dalcin ierr = VecGetArray(b,&bv);CHKERRQ(ierr); 159*16d9e3a6SLisandro Dalcin ierr = VecGetArray(x,&xv);CHKERRQ(ierr); 160*16d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,bv,sbv); 161*16d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,xv,sxv); 162*16d9e3a6SLisandro Dalcin 163*16d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 164*16d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr); 165*16d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr); 166*16d9e3a6SLisandro Dalcin hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv); 167*16d9e3a6SLisandro Dalcin /* error code of 1 in BoomerAMG merely means convergence not achieved */ 168*16d9e3a6SLisandro Dalcin if (hierr && (hierr != 1 || jac->solve != HYPRE_BoomerAMGSolve)) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 169*16d9e3a6SLisandro Dalcin if (hierr) hypre__global_error = 0; 170*16d9e3a6SLisandro Dalcin 171*16d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,sbv,bv); 172*16d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,sxv,xv); 173*16d9e3a6SLisandro Dalcin ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr); 174*16d9e3a6SLisandro Dalcin ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr); 175*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 176*16d9e3a6SLisandro Dalcin } 177*16d9e3a6SLisandro Dalcin 178*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 179*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCDestroy_HYPRE" 180*16d9e3a6SLisandro Dalcin static PetscErrorCode PCDestroy_HYPRE(PC pc) 181*16d9e3a6SLisandro Dalcin { 182*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 183*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 184*16d9e3a6SLisandro Dalcin 185*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 186*16d9e3a6SLisandro Dalcin if (jac->ij) { ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); } 187*16d9e3a6SLisandro Dalcin if (jac->b) { ierr = HYPRE_IJVectorDestroy(jac->b);CHKERRQ(ierr); } 188*16d9e3a6SLisandro Dalcin if (jac->x) { ierr = HYPRE_IJVectorDestroy(jac->x);CHKERRQ(ierr); } 189*16d9e3a6SLisandro Dalcin if (jac->destroy) { ierr = (*jac->destroy)(jac->hsolver);CHKERRQ(ierr); } 190*16d9e3a6SLisandro Dalcin ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr); 191*16d9e3a6SLisandro Dalcin if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);} 192*16d9e3a6SLisandro Dalcin ierr = PetscFree(jac);CHKERRQ(ierr); 193*16d9e3a6SLisandro Dalcin 194*16d9e3a6SLisandro Dalcin ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr); 195*16d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);CHKERRQ(ierr); 196*16d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);CHKERRQ(ierr); 197*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 198*16d9e3a6SLisandro Dalcin } 199*16d9e3a6SLisandro Dalcin 200*16d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 201*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 202*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut" 203*16d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc) 204*16d9e3a6SLisandro Dalcin { 205*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 206*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 207*16d9e3a6SLisandro Dalcin PetscTruth flag; 208*16d9e3a6SLisandro Dalcin 209*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 210*16d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr); 211*16d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr); 212*16d9e3a6SLisandro Dalcin if (flag) { 213*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 214*16d9e3a6SLisandro Dalcin } 215*16d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr); 216*16d9e3a6SLisandro Dalcin if (flag) { 217*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutSetDropTolerance(jac->hsolver,jac->tol);CHKERRQ(ierr); 218*16d9e3a6SLisandro Dalcin } 219*16d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr); 220*16d9e3a6SLisandro Dalcin if (flag) { 221*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutSetFactorRowSize(jac->hsolver,jac->factorrowsize);CHKERRQ(ierr); 222*16d9e3a6SLisandro Dalcin } 223*16d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 224*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 225*16d9e3a6SLisandro Dalcin } 226*16d9e3a6SLisandro Dalcin 227*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 228*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Pilut" 229*16d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer) 230*16d9e3a6SLisandro Dalcin { 231*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 232*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 233*16d9e3a6SLisandro Dalcin PetscTruth iascii; 234*16d9e3a6SLisandro Dalcin 235*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 236*16d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 237*16d9e3a6SLisandro Dalcin if (iascii) { 238*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut preconditioning\n");CHKERRQ(ierr); 239*16d9e3a6SLisandro Dalcin if (jac->maxiter != PETSC_DEFAULT) { 240*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr); 241*16d9e3a6SLisandro Dalcin } else { 242*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr); 243*16d9e3a6SLisandro Dalcin } 244*16d9e3a6SLisandro Dalcin if (jac->tol != PETSC_DEFAULT) { 245*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr); 246*16d9e3a6SLisandro Dalcin } else { 247*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr); 248*16d9e3a6SLisandro Dalcin } 249*16d9e3a6SLisandro Dalcin if (jac->factorrowsize != PETSC_DEFAULT) { 250*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr); 251*16d9e3a6SLisandro Dalcin } else { 252*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default factor row size \n");CHKERRQ(ierr); 253*16d9e3a6SLisandro Dalcin } 254*16d9e3a6SLisandro Dalcin } 255*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 256*16d9e3a6SLisandro Dalcin } 257*16d9e3a6SLisandro Dalcin 258*16d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 259*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 260*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid" 261*16d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc) 262*16d9e3a6SLisandro Dalcin { 263*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 264*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 265*16d9e3a6SLisandro Dalcin PetscTruth flag; 266*16d9e3a6SLisandro Dalcin char *args[2]; 267*16d9e3a6SLisandro Dalcin 268*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 269*16d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr); 270*16d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr); 271*16d9e3a6SLisandro Dalcin if (flag) { 272*16d9e3a6SLisandro Dalcin char levels[16]; 273*16d9e3a6SLisandro Dalcin if (jac->levels < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels); 274*16d9e3a6SLisandro Dalcin sprintf(levels,"%d",jac->levels); 275*16d9e3a6SLisandro Dalcin args[0] = (char*)"-level"; args[1] = levels; 276*16d9e3a6SLisandro Dalcin ierr = HYPRE_EuclidSetParams(jac->hsolver,2,args);CHKERRQ(ierr); 277*16d9e3a6SLisandro Dalcin } 278*16d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);CHKERRQ(ierr); 279*16d9e3a6SLisandro Dalcin if (jac->bjilu) { 280*16d9e3a6SLisandro Dalcin args[0] =(char*) "-bj"; args[1] = (char*)"1"; 281*16d9e3a6SLisandro Dalcin ierr = HYPRE_EuclidSetParams(jac->hsolver,2,args);CHKERRQ(ierr); 282*16d9e3a6SLisandro Dalcin } 283*16d9e3a6SLisandro Dalcin 284*16d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);CHKERRQ(ierr); 285*16d9e3a6SLisandro Dalcin if (jac->printstatistics) { 286*16d9e3a6SLisandro Dalcin args[0] = (char*)"-eu_stats"; args[1] = (char*)"1"; 287*16d9e3a6SLisandro Dalcin ierr = HYPRE_EuclidSetParams(jac->hsolver,2,args);CHKERRQ(ierr); 288*16d9e3a6SLisandro Dalcin args[0] = (char*)"-eu_mem"; args[1] = (char*)"1"; 289*16d9e3a6SLisandro Dalcin ierr = HYPRE_EuclidSetParams(jac->hsolver,2,args);CHKERRQ(ierr); 290*16d9e3a6SLisandro Dalcin } 291*16d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 292*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 293*16d9e3a6SLisandro Dalcin } 294*16d9e3a6SLisandro Dalcin 295*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 296*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Euclid" 297*16d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer) 298*16d9e3a6SLisandro Dalcin { 299*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 300*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 301*16d9e3a6SLisandro Dalcin PetscTruth iascii; 302*16d9e3a6SLisandro Dalcin 303*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 304*16d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 305*16d9e3a6SLisandro Dalcin if (iascii) { 306*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid preconditioning\n");CHKERRQ(ierr); 307*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr); 308*16d9e3a6SLisandro Dalcin if (jac->bjilu) { 309*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr); 310*16d9e3a6SLisandro Dalcin } 311*16d9e3a6SLisandro Dalcin } 312*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 313*16d9e3a6SLisandro Dalcin } 314*16d9e3a6SLisandro Dalcin 315*16d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 316*16d9e3a6SLisandro Dalcin 317*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 318*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG" 319*16d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x) 320*16d9e3a6SLisandro Dalcin { 321*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 322*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 323*16d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 324*16d9e3a6SLisandro Dalcin PetscScalar *bv,*xv; 325*16d9e3a6SLisandro Dalcin HYPRE_ParVector jbv,jxv; 326*16d9e3a6SLisandro Dalcin PetscScalar *sbv,*sxv; 327*16d9e3a6SLisandro Dalcin int hierr; 328*16d9e3a6SLisandro Dalcin 329*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 330*16d9e3a6SLisandro Dalcin ierr = VecSet(x,0.0);CHKERRQ(ierr); 331*16d9e3a6SLisandro Dalcin ierr = VecGetArray(b,&bv);CHKERRQ(ierr); 332*16d9e3a6SLisandro Dalcin ierr = VecGetArray(x,&xv);CHKERRQ(ierr); 333*16d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,bv,sbv); 334*16d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,xv,sxv); 335*16d9e3a6SLisandro Dalcin 336*16d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 337*16d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr); 338*16d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr); 339*16d9e3a6SLisandro Dalcin 340*16d9e3a6SLisandro Dalcin hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv); 341*16d9e3a6SLisandro Dalcin /* error code of 1 in BoomerAMG merely means convergence not achieved */ 342*16d9e3a6SLisandro Dalcin if (hierr && (hierr != 1)) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 343*16d9e3a6SLisandro Dalcin if (hierr) hypre__global_error = 0; 344*16d9e3a6SLisandro Dalcin 345*16d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,sbv,bv); 346*16d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,sxv,xv); 347*16d9e3a6SLisandro Dalcin ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr); 348*16d9e3a6SLisandro Dalcin ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr); 349*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 350*16d9e3a6SLisandro Dalcin } 351*16d9e3a6SLisandro Dalcin 352*16d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCycleType[] = {"","V","W"}; 353*16d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout"}; 354*16d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGMeasureType[] = {"local","global"}; 355*16d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGRelaxType[] = {"Jacobi","sequential-Gauss-Seidel","","SOR/Jacobi","backward-SOR/Jacobi","","symmetric-SOR/Jacobi","","","Gaussian-elimination"}; 356*16d9e3a6SLisandro Dalcin 357*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 358*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG" 359*16d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc) 360*16d9e3a6SLisandro Dalcin { 361*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 362*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 363*16d9e3a6SLisandro Dalcin int n,indx; 364*16d9e3a6SLisandro Dalcin PetscTruth flg, tmp_truth; 365*16d9e3a6SLisandro Dalcin double tmpdbl, twodbl[2]; 366*16d9e3a6SLisandro Dalcin 367*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 368*16d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr); 369*16d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr); 370*16d9e3a6SLisandro Dalcin if (flg) { 371*16d9e3a6SLisandro Dalcin jac->cycletype = indx; 372*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr); 373*16d9e3a6SLisandro Dalcin } 374*16d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr); 375*16d9e3a6SLisandro Dalcin if (flg) { 376*16d9e3a6SLisandro Dalcin if (jac->maxlevels < 2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels); 377*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr); 378*16d9e3a6SLisandro Dalcin } 379*16d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr); 380*16d9e3a6SLisandro Dalcin if (flg) { 381*16d9e3a6SLisandro Dalcin if (jac->maxiter < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter); 382*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 383*16d9e3a6SLisandro Dalcin } 384*16d9e3a6SLisandro Dalcin ierr = PetscOptionsScalar("-pc_hypre_boomeramg_tol","Convergence tolerance PER hypre call","None",jac->tol,&jac->tol,&flg);CHKERRQ(ierr); 385*16d9e3a6SLisandro Dalcin if (flg) { 386*16d9e3a6SLisandro Dalcin if (jac->tol < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be great than or equal zero",jac->tol); 387*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 388*16d9e3a6SLisandro Dalcin } 389*16d9e3a6SLisandro Dalcin 390*16d9e3a6SLisandro Dalcin ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr); 391*16d9e3a6SLisandro Dalcin if (flg) { 392*16d9e3a6SLisandro Dalcin if (jac->truncfactor < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor); 393*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr); 394*16d9e3a6SLisandro Dalcin } 395*16d9e3a6SLisandro Dalcin 396*16d9e3a6SLisandro Dalcin ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr); 397*16d9e3a6SLisandro Dalcin if (flg) { 398*16d9e3a6SLisandro Dalcin if (jac->strongthreshold < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold); 399*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr); 400*16d9e3a6SLisandro Dalcin } 401*16d9e3a6SLisandro Dalcin ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr); 402*16d9e3a6SLisandro Dalcin if (flg) { 403*16d9e3a6SLisandro Dalcin if (jac->maxrowsum < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum); 404*16d9e3a6SLisandro Dalcin if (jac->maxrowsum > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum); 405*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr); 406*16d9e3a6SLisandro Dalcin } 407*16d9e3a6SLisandro Dalcin 408*16d9e3a6SLisandro Dalcin /* Grid sweeps */ 409*16d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_all","Number of sweeps for all grid levels (fine, up, and down)","None",jac->gridsweeps[0],&indx,&flg);CHKERRQ(ierr); 410*16d9e3a6SLisandro Dalcin if (flg) { 411*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver,indx);CHKERRQ(ierr); 412*16d9e3a6SLisandro Dalcin /* modify the jac structure so we can view the updated options with PC_View */ 413*16d9e3a6SLisandro Dalcin jac->gridsweeps[0] = indx; 414*16d9e3a6SLisandro Dalcin jac->gridsweeps[1] = jac->gridsweeps[2] = jac->gridsweeps[0]; 415*16d9e3a6SLisandro Dalcin jac->gridsweeps[3] = 1; /*The coarse level is not affected by this function - hypre code sets to 1*/ 416*16d9e3a6SLisandro Dalcin } 417*16d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_fine","Number of sweeps for the fine level","None",jac->gridsweeps[0],&indx,&flg);CHKERRQ(ierr); 418*16d9e3a6SLisandro Dalcin if (flg) { 419*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 0);CHKERRQ(ierr); 420*16d9e3a6SLisandro Dalcin jac->gridsweeps[0] = indx; 421*16d9e3a6SLisandro Dalcin } 422*16d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr); 423*16d9e3a6SLisandro Dalcin if (flg) { 424*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 1);CHKERRQ(ierr); 425*16d9e3a6SLisandro Dalcin jac->gridsweeps[1] = indx; 426*16d9e3a6SLisandro Dalcin } 427*16d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr); 428*16d9e3a6SLisandro Dalcin if (flg) { 429*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 2);CHKERRQ(ierr); 430*16d9e3a6SLisandro Dalcin jac->gridsweeps[2] = indx; 431*16d9e3a6SLisandro Dalcin } 432*16d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[3],&indx,&flg);CHKERRQ(ierr); 433*16d9e3a6SLisandro Dalcin if (flg) { 434*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 3);CHKERRQ(ierr); 435*16d9e3a6SLisandro Dalcin jac->gridsweeps[3] = indx; 436*16d9e3a6SLisandro Dalcin } 437*16d9e3a6SLisandro Dalcin 438*16d9e3a6SLisandro Dalcin /* Relax type */ 439*16d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_all","Relax type for fine, up, and down cycles (coarse level set to gaussian elimination)","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);CHKERRQ(ierr); 440*16d9e3a6SLisandro Dalcin if (flg) { 441*16d9e3a6SLisandro Dalcin jac->relaxtype[0] = jac->relaxtype[1] = jac->relaxtype[2] = indx; 442*16d9e3a6SLisandro Dalcin jac->relaxtype[3] = 9; /* hypre code sets coarse grid to 9 (G.E.)*/ 443*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, indx);CHKERRQ(ierr); 444*16d9e3a6SLisandro Dalcin } 445*16d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_fine","Relax type on fine grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);CHKERRQ(ierr); 446*16d9e3a6SLisandro Dalcin if (flg) { 447*16d9e3a6SLisandro Dalcin jac->relaxtype[0] = indx; 448*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 0);CHKERRQ(ierr); 449*16d9e3a6SLisandro Dalcin } 450*16d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);CHKERRQ(ierr); 451*16d9e3a6SLisandro Dalcin if (flg) { 452*16d9e3a6SLisandro Dalcin jac->relaxtype[1] = indx; 453*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 1);CHKERRQ(ierr); 454*16d9e3a6SLisandro Dalcin } 455*16d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[3],&indx,&flg);CHKERRQ(ierr); 456*16d9e3a6SLisandro Dalcin if (flg) { 457*16d9e3a6SLisandro Dalcin jac->relaxtype[2] = indx; 458*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 2);CHKERRQ(ierr); 459*16d9e3a6SLisandro Dalcin } 460*16d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr); 461*16d9e3a6SLisandro Dalcin if (flg) { 462*16d9e3a6SLisandro Dalcin jac->relaxtype[3] = indx; 463*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 3);CHKERRQ(ierr); 464*16d9e3a6SLisandro Dalcin } 465*16d9e3a6SLisandro Dalcin 466*16d9e3a6SLisandro Dalcin /* Relaxation Weight */ 467*16d9e3a6SLisandro 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); 468*16d9e3a6SLisandro Dalcin if (flg) { 469*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxWt(jac->hsolver,tmpdbl);CHKERRQ(ierr); 470*16d9e3a6SLisandro Dalcin jac->relaxweight = tmpdbl; 471*16d9e3a6SLisandro Dalcin } 472*16d9e3a6SLisandro Dalcin 473*16d9e3a6SLisandro Dalcin n=2; 474*16d9e3a6SLisandro Dalcin twodbl[0] = twodbl[1] = 1.0; 475*16d9e3a6SLisandro 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); 476*16d9e3a6SLisandro Dalcin if (flg) { 477*16d9e3a6SLisandro Dalcin if (n == 2) { 478*16d9e3a6SLisandro Dalcin indx = (int)PetscAbsReal(twodbl[1]); 479*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetLevelRelaxWt(jac->hsolver,twodbl[0],indx);CHKERRQ(ierr); 480*16d9e3a6SLisandro Dalcin } else { 481*16d9e3a6SLisandro Dalcin SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Relax weight level: you must provide 2 values separated by a comma (and no space), you provided %d",n); 482*16d9e3a6SLisandro Dalcin } 483*16d9e3a6SLisandro Dalcin } 484*16d9e3a6SLisandro Dalcin 485*16d9e3a6SLisandro Dalcin /* Outer relaxation Weight */ 486*16d9e3a6SLisandro 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); 487*16d9e3a6SLisandro Dalcin if (flg) { 488*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetOuterWt( jac->hsolver, tmpdbl);CHKERRQ(ierr); 489*16d9e3a6SLisandro Dalcin jac->outerrelaxweight = tmpdbl; 490*16d9e3a6SLisandro Dalcin } 491*16d9e3a6SLisandro Dalcin 492*16d9e3a6SLisandro Dalcin n=2; 493*16d9e3a6SLisandro Dalcin twodbl[0] = twodbl[1] = 1.0; 494*16d9e3a6SLisandro 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); 495*16d9e3a6SLisandro Dalcin if (flg) { 496*16d9e3a6SLisandro Dalcin if (n == 2) { 497*16d9e3a6SLisandro Dalcin indx = (int)PetscAbsReal(twodbl[1]); 498*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetLevelOuterWt( jac->hsolver, twodbl[0], indx);CHKERRQ(ierr); 499*16d9e3a6SLisandro Dalcin } else { 500*16d9e3a6SLisandro Dalcin SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Relax weight outer level: You must provide 2 values separated by a comma (and no space), you provided %d",n); 501*16d9e3a6SLisandro Dalcin } 502*16d9e3a6SLisandro Dalcin } 503*16d9e3a6SLisandro Dalcin 504*16d9e3a6SLisandro Dalcin /* the Relax Order */ 505*16d9e3a6SLisandro Dalcin /* ierr = PetscOptionsName("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", &flg);CHKERRQ(ierr); */ 506*16d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr); 507*16d9e3a6SLisandro Dalcin 508*16d9e3a6SLisandro Dalcin if (flg) { 509*16d9e3a6SLisandro Dalcin jac->relaxorder = 0; 510*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr); 511*16d9e3a6SLisandro Dalcin } 512*16d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr); 513*16d9e3a6SLisandro Dalcin if (flg) { 514*16d9e3a6SLisandro Dalcin jac->measuretype = indx; 515*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr); 516*16d9e3a6SLisandro Dalcin } 517*16d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,7,HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr); 518*16d9e3a6SLisandro Dalcin if (flg) { 519*16d9e3a6SLisandro Dalcin jac->coarsentype = indx; 520*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr); 521*16d9e3a6SLisandro Dalcin } 522*16d9e3a6SLisandro Dalcin ierr = PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);CHKERRQ(ierr); 523*16d9e3a6SLisandro Dalcin if (flg) { 524*16d9e3a6SLisandro Dalcin int level=3; 525*16d9e3a6SLisandro Dalcin jac->printstatistics = PETSC_TRUE; 526*16d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);CHKERRQ(ierr); 527*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetPrintLevel(jac->hsolver,level);CHKERRQ(ierr); 528*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetDebugFlag(jac->hsolver,level);CHKERRQ(ierr); 529*16d9e3a6SLisandro Dalcin } 530*16d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 531*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 532*16d9e3a6SLisandro Dalcin } 533*16d9e3a6SLisandro Dalcin 534*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 535*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG" 536*16d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its) 537*16d9e3a6SLisandro Dalcin { 538*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 539*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 540*16d9e3a6SLisandro Dalcin 541*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 542*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,its*jac->maxiter);CHKERRQ(ierr); 543*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,PetscMin(rtol,jac->tol));CHKERRQ(ierr); 544*16d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_TRUE; 545*16d9e3a6SLisandro Dalcin ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr); 546*16d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_FALSE; 547*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 548*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 549*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 550*16d9e3a6SLisandro Dalcin } 551*16d9e3a6SLisandro Dalcin 552*16d9e3a6SLisandro Dalcin 553*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 554*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_BoomerAMG" 555*16d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer) 556*16d9e3a6SLisandro Dalcin { 557*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 558*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 559*16d9e3a6SLisandro Dalcin PetscTruth iascii; 560*16d9e3a6SLisandro Dalcin 561*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 562*16d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 563*16d9e3a6SLisandro Dalcin if (iascii) { 564*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr); 565*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr); 566*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr); 567*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr); 568*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr); 569*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr); 570*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr); 571*16d9e3a6SLisandro Dalcin 572*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps on fine grid %d\n",jac->gridsweeps[0]);CHKERRQ(ierr); 573*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps down %d\n",jac->gridsweeps[1]);CHKERRQ(ierr); 574*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps up %d\n",jac->gridsweeps[2]);CHKERRQ(ierr); 575*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps on coarse %d\n",jac->gridsweeps[3]);CHKERRQ(ierr); 576*16d9e3a6SLisandro Dalcin 577*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax on fine grid %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr); 578*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax down %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr); 579*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax up %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr); 580*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax on coarse %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[3]]);CHKERRQ(ierr); 581*16d9e3a6SLisandro Dalcin 582*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax weight (all) %G\n",jac->relaxweight);CHKERRQ(ierr); 583*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr); 584*16d9e3a6SLisandro Dalcin 585*16d9e3a6SLisandro Dalcin if (jac->relaxorder) { 586*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr); 587*16d9e3a6SLisandro Dalcin } else { 588*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr); 589*16d9e3a6SLisandro Dalcin } 590*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Measure type %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr); 591*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Coarsen type %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr); 592*16d9e3a6SLisandro Dalcin } 593*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 594*16d9e3a6SLisandro Dalcin } 595*16d9e3a6SLisandro Dalcin 596*16d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 597*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 598*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails" 599*16d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc) 600*16d9e3a6SLisandro Dalcin { 601*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 602*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 603*16d9e3a6SLisandro Dalcin int indx; 604*16d9e3a6SLisandro Dalcin PetscTruth flag; 605*16d9e3a6SLisandro Dalcin const char *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"}; 606*16d9e3a6SLisandro Dalcin 607*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 608*16d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr); 609*16d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr); 610*16d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr); 611*16d9e3a6SLisandro Dalcin if (flag) { 612*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr); 613*16d9e3a6SLisandro Dalcin } 614*16d9e3a6SLisandro Dalcin 615*16d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr); 616*16d9e3a6SLisandro Dalcin if (flag) { 617*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr); 618*16d9e3a6SLisandro Dalcin } 619*16d9e3a6SLisandro Dalcin 620*16d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr); 621*16d9e3a6SLisandro Dalcin if (flag) { 622*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr); 623*16d9e3a6SLisandro Dalcin } 624*16d9e3a6SLisandro Dalcin 625*16d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,&flag);CHKERRQ(ierr); 626*16d9e3a6SLisandro Dalcin if (flag) { 627*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr); 628*16d9e3a6SLisandro Dalcin } 629*16d9e3a6SLisandro Dalcin 630*16d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,&flag);CHKERRQ(ierr); 631*16d9e3a6SLisandro Dalcin if (flag) { 632*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr); 633*16d9e3a6SLisandro Dalcin } 634*16d9e3a6SLisandro Dalcin 635*16d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);CHKERRQ(ierr); 636*16d9e3a6SLisandro Dalcin if (flag) { 637*16d9e3a6SLisandro Dalcin jac->symt = indx; 638*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr); 639*16d9e3a6SLisandro Dalcin } 640*16d9e3a6SLisandro Dalcin 641*16d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 642*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 643*16d9e3a6SLisandro Dalcin } 644*16d9e3a6SLisandro Dalcin 645*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 646*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_ParaSails" 647*16d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer) 648*16d9e3a6SLisandro Dalcin { 649*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 650*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 651*16d9e3a6SLisandro Dalcin PetscTruth iascii; 652*16d9e3a6SLisandro Dalcin const char *symt = 0;; 653*16d9e3a6SLisandro Dalcin 654*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 655*16d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 656*16d9e3a6SLisandro Dalcin if (iascii) { 657*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails preconditioning\n");CHKERRQ(ierr); 658*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr); 659*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr); 660*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr); 661*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr); 662*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: reuse nonzero structure %s\n",PetscTruths[jac->ruse]);CHKERRQ(ierr); 663*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: print info to screen %s\n",PetscTruths[jac->logging]);CHKERRQ(ierr); 664*16d9e3a6SLisandro Dalcin if (!jac->symt) { 665*16d9e3a6SLisandro Dalcin symt = "nonsymmetric matrix and preconditioner"; 666*16d9e3a6SLisandro Dalcin } else if (jac->symt == 1) { 667*16d9e3a6SLisandro Dalcin symt = "SPD matrix and preconditioner"; 668*16d9e3a6SLisandro Dalcin } else if (jac->symt == 2) { 669*16d9e3a6SLisandro Dalcin symt = "nonsymmetric matrix but SPD preconditioner"; 670*16d9e3a6SLisandro Dalcin } else { 671*16d9e3a6SLisandro Dalcin SETERRQ1(PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt); 672*16d9e3a6SLisandro Dalcin } 673*16d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr); 674*16d9e3a6SLisandro Dalcin } 675*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 676*16d9e3a6SLisandro Dalcin } 677*16d9e3a6SLisandro Dalcin /* ---------------------------------------------------------------------------------*/ 678*16d9e3a6SLisandro Dalcin 679*16d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 680*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 681*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType_HYPRE" 682*16d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType_HYPRE(PC pc,const char *name[]) 683*16d9e3a6SLisandro Dalcin { 684*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 685*16d9e3a6SLisandro Dalcin 686*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 687*16d9e3a6SLisandro Dalcin *name = jac->hypre_type; 688*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 689*16d9e3a6SLisandro Dalcin } 690*16d9e3a6SLisandro Dalcin EXTERN_C_END 691*16d9e3a6SLisandro Dalcin 692*16d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 693*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 694*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType_HYPRE" 695*16d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType_HYPRE(PC pc,const char name[]) 696*16d9e3a6SLisandro Dalcin { 697*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 698*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 699*16d9e3a6SLisandro Dalcin PetscTruth flag; 700*16d9e3a6SLisandro Dalcin 701*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 702*16d9e3a6SLisandro Dalcin if (jac->hypre_type) { 703*16d9e3a6SLisandro Dalcin ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr); 704*16d9e3a6SLisandro Dalcin if (!flag) SETERRQ(PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set"); 705*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 706*16d9e3a6SLisandro Dalcin } else { 707*16d9e3a6SLisandro Dalcin ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr); 708*16d9e3a6SLisandro Dalcin } 709*16d9e3a6SLisandro Dalcin 710*16d9e3a6SLisandro Dalcin jac->maxiter = PETSC_DEFAULT; 711*16d9e3a6SLisandro Dalcin jac->tol = PETSC_DEFAULT; 712*16d9e3a6SLisandro Dalcin jac->printstatistics = PetscLogPrintInfo; 713*16d9e3a6SLisandro Dalcin 714*16d9e3a6SLisandro Dalcin ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr); 715*16d9e3a6SLisandro Dalcin if (flag) { 716*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutCreate(jac->comm_hypre,&jac->hsolver); 717*16d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut; 718*16d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_Pilut; 719*16d9e3a6SLisandro Dalcin jac->destroy = HYPRE_ParCSRPilutDestroy; 720*16d9e3a6SLisandro Dalcin jac->setup = HYPRE_ParCSRPilutSetup; 721*16d9e3a6SLisandro Dalcin jac->solve = HYPRE_ParCSRPilutSolve; 722*16d9e3a6SLisandro Dalcin jac->factorrowsize = PETSC_DEFAULT; 723*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 724*16d9e3a6SLisandro Dalcin } 725*16d9e3a6SLisandro Dalcin ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr); 726*16d9e3a6SLisandro Dalcin if (flag) { 727*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsCreate(jac->comm_hypre,&jac->hsolver); 728*16d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails; 729*16d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_ParaSails; 730*16d9e3a6SLisandro Dalcin jac->destroy = HYPRE_ParaSailsDestroy; 731*16d9e3a6SLisandro Dalcin jac->setup = HYPRE_ParaSailsSetup; 732*16d9e3a6SLisandro Dalcin jac->solve = HYPRE_ParaSailsSolve; 733*16d9e3a6SLisandro Dalcin /* initialize */ 734*16d9e3a6SLisandro Dalcin jac->nlevels = 1; 735*16d9e3a6SLisandro Dalcin jac->threshhold = .1; 736*16d9e3a6SLisandro Dalcin jac->filter = .1; 737*16d9e3a6SLisandro Dalcin jac->loadbal = 0; 738*16d9e3a6SLisandro Dalcin if (PetscLogPrintInfo) { 739*16d9e3a6SLisandro Dalcin jac->logging = (int) PETSC_TRUE; 740*16d9e3a6SLisandro Dalcin } else { 741*16d9e3a6SLisandro Dalcin jac->logging = (int) PETSC_FALSE; 742*16d9e3a6SLisandro Dalcin } 743*16d9e3a6SLisandro Dalcin jac->ruse = (int) PETSC_FALSE; 744*16d9e3a6SLisandro Dalcin jac->symt = 0; 745*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr); 746*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr); 747*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr); 748*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr); 749*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr); 750*16d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr); 751*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 752*16d9e3a6SLisandro Dalcin } 753*16d9e3a6SLisandro Dalcin ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr); 754*16d9e3a6SLisandro Dalcin if (flag) { 755*16d9e3a6SLisandro Dalcin ierr = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver); 756*16d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid; 757*16d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_Euclid; 758*16d9e3a6SLisandro Dalcin jac->destroy = HYPRE_EuclidDestroy; 759*16d9e3a6SLisandro Dalcin jac->setup = HYPRE_EuclidSetup; 760*16d9e3a6SLisandro Dalcin jac->solve = HYPRE_EuclidSolve; 761*16d9e3a6SLisandro Dalcin /* initialization */ 762*16d9e3a6SLisandro Dalcin jac->bjilu = PETSC_FALSE; 763*16d9e3a6SLisandro Dalcin jac->levels = 1; 764*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 765*16d9e3a6SLisandro Dalcin } 766*16d9e3a6SLisandro Dalcin ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr); 767*16d9e3a6SLisandro Dalcin if (flag) { 768*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGCreate(&jac->hsolver); 769*16d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_BoomerAMG; 770*16d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_BoomerAMG; 771*16d9e3a6SLisandro Dalcin pc->ops->applytranspose = PCApplyTranspose_HYPRE_BoomerAMG; 772*16d9e3a6SLisandro Dalcin pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG; 773*16d9e3a6SLisandro Dalcin jac->destroy = HYPRE_BoomerAMGDestroy; 774*16d9e3a6SLisandro Dalcin jac->setup = HYPRE_BoomerAMGSetup; 775*16d9e3a6SLisandro Dalcin jac->solve = HYPRE_BoomerAMGSolve; 776*16d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_FALSE; 777*16d9e3a6SLisandro Dalcin /* these defaults match the hypre defaults */ 778*16d9e3a6SLisandro Dalcin jac->cycletype = 1; 779*16d9e3a6SLisandro Dalcin jac->maxlevels = 25; 780*16d9e3a6SLisandro Dalcin jac->maxiter = 1; 781*16d9e3a6SLisandro Dalcin jac->tol = 1.e-7; 782*16d9e3a6SLisandro Dalcin jac->truncfactor = 0.0; 783*16d9e3a6SLisandro Dalcin jac->strongthreshold = .25; 784*16d9e3a6SLisandro Dalcin jac->maxrowsum = .9; 785*16d9e3a6SLisandro Dalcin jac->coarsentype = 6; 786*16d9e3a6SLisandro Dalcin jac->measuretype = 0; 787*16d9e3a6SLisandro Dalcin jac->gridsweeps[0] = jac->gridsweeps[1] = jac->gridsweeps[2] = jac->gridsweeps[3] = 1; 788*16d9e3a6SLisandro Dalcin jac->relaxtype[0] = jac->relaxtype[1] = jac->relaxtype[2] = 3; 789*16d9e3a6SLisandro Dalcin jac->relaxtype[3] = 9; 790*16d9e3a6SLisandro Dalcin jac->relaxweight = 1.0; 791*16d9e3a6SLisandro Dalcin jac->outerrelaxweight = 1.0; 792*16d9e3a6SLisandro Dalcin jac->relaxorder = 1; 793*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr); 794*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr); 795*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 796*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 797*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr); 798*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr); 799*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr); 800*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr); 801*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr); 802*16d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr); 803*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 804*16d9e3a6SLisandro Dalcin } 805*16d9e3a6SLisandro Dalcin ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr); 806*16d9e3a6SLisandro Dalcin jac->hypre_type = PETSC_NULL; 807*16d9e3a6SLisandro Dalcin SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name); 808*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 809*16d9e3a6SLisandro Dalcin } 810*16d9e3a6SLisandro Dalcin EXTERN_C_END 811*16d9e3a6SLisandro Dalcin 812*16d9e3a6SLisandro Dalcin /* 813*16d9e3a6SLisandro Dalcin It only gets here if the HYPRE type has not been set before the call to 814*16d9e3a6SLisandro Dalcin ...SetFromOptions() which actually is most of the time 815*16d9e3a6SLisandro Dalcin */ 816*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 817*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE" 818*16d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE(PC pc) 819*16d9e3a6SLisandro Dalcin { 820*16d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 821*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 822*16d9e3a6SLisandro Dalcin int indx; 823*16d9e3a6SLisandro Dalcin const char *type[] = {"pilut","parasails","boomeramg","euclid"}; 824*16d9e3a6SLisandro Dalcin PetscTruth flg; 825*16d9e3a6SLisandro Dalcin 826*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 827*16d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr); 828*16d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"pilut",&indx,&flg);CHKERRQ(ierr); 829*16d9e3a6SLisandro Dalcin if (PetscOptionsPublishCount) { /* force the default if it was not yet set and user did not set with option */ 830*16d9e3a6SLisandro Dalcin if (!flg && !jac->hypre_type) { 831*16d9e3a6SLisandro Dalcin flg = PETSC_TRUE; 832*16d9e3a6SLisandro Dalcin indx = 0; 833*16d9e3a6SLisandro Dalcin } 834*16d9e3a6SLisandro Dalcin } 835*16d9e3a6SLisandro Dalcin if (flg) { 836*16d9e3a6SLisandro Dalcin ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr); 837*16d9e3a6SLisandro Dalcin } 838*16d9e3a6SLisandro Dalcin if (pc->ops->setfromoptions) { 839*16d9e3a6SLisandro Dalcin ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr); 840*16d9e3a6SLisandro Dalcin } 841*16d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 842*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 843*16d9e3a6SLisandro Dalcin } 844*16d9e3a6SLisandro Dalcin 845*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 846*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType" 847*16d9e3a6SLisandro Dalcin /*@C 848*16d9e3a6SLisandro Dalcin PCHYPRESetType - Sets which hypre preconditioner you wish to use 849*16d9e3a6SLisandro Dalcin 850*16d9e3a6SLisandro Dalcin Input Parameters: 851*16d9e3a6SLisandro Dalcin + pc - the preconditioner context 852*16d9e3a6SLisandro Dalcin - name - either pilut, parasails, boomeramg, euclid 853*16d9e3a6SLisandro Dalcin 854*16d9e3a6SLisandro Dalcin Options Database Keys: 855*16d9e3a6SLisandro Dalcin -pc_hypre_type - One of pilut, parasails, boomeramg, euclid 856*16d9e3a6SLisandro Dalcin 857*16d9e3a6SLisandro Dalcin Level: intermediate 858*16d9e3a6SLisandro Dalcin 859*16d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 860*16d9e3a6SLisandro Dalcin PCHYPRE 861*16d9e3a6SLisandro Dalcin 862*16d9e3a6SLisandro Dalcin @*/ 863*16d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC pc,const char name[]) 864*16d9e3a6SLisandro Dalcin { 865*16d9e3a6SLisandro Dalcin PetscErrorCode ierr,(*f)(PC,const char[]); 866*16d9e3a6SLisandro Dalcin 867*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 868*16d9e3a6SLisandro Dalcin PetscValidHeaderSpecific(pc,PC_COOKIE,1); 869*16d9e3a6SLisandro Dalcin PetscValidCharPointer(name,2); 870*16d9e3a6SLisandro Dalcin ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPRESetType_C",(void (**)(void))&f);CHKERRQ(ierr); 871*16d9e3a6SLisandro Dalcin if (f) { 872*16d9e3a6SLisandro Dalcin ierr = (*f)(pc,name);CHKERRQ(ierr); 873*16d9e3a6SLisandro Dalcin } 874*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 875*16d9e3a6SLisandro Dalcin } 876*16d9e3a6SLisandro Dalcin 877*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 878*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType" 879*16d9e3a6SLisandro Dalcin /*@C 880*16d9e3a6SLisandro Dalcin PCHYPREGetType - Gets which hypre preconditioner you are using 881*16d9e3a6SLisandro Dalcin 882*16d9e3a6SLisandro Dalcin Input Parameter: 883*16d9e3a6SLisandro Dalcin . pc - the preconditioner context 884*16d9e3a6SLisandro Dalcin 885*16d9e3a6SLisandro Dalcin Output Parameter: 886*16d9e3a6SLisandro Dalcin . name - either pilut, parasails, boomeramg, euclid 887*16d9e3a6SLisandro Dalcin 888*16d9e3a6SLisandro Dalcin Level: intermediate 889*16d9e3a6SLisandro Dalcin 890*16d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC, 891*16d9e3a6SLisandro Dalcin PCHYPRE 892*16d9e3a6SLisandro Dalcin 893*16d9e3a6SLisandro Dalcin @*/ 894*16d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC pc,const char *name[]) 895*16d9e3a6SLisandro Dalcin { 896*16d9e3a6SLisandro Dalcin PetscErrorCode ierr,(*f)(PC,const char*[]); 897*16d9e3a6SLisandro Dalcin 898*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 899*16d9e3a6SLisandro Dalcin PetscValidHeaderSpecific(pc,PC_COOKIE,1); 900*16d9e3a6SLisandro Dalcin PetscValidPointer(name,2); 901*16d9e3a6SLisandro Dalcin ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPREGetType_C",(void (**)(void))&f);CHKERRQ(ierr); 902*16d9e3a6SLisandro Dalcin if (f) { 903*16d9e3a6SLisandro Dalcin ierr = (*f)(pc,name);CHKERRQ(ierr); 904*16d9e3a6SLisandro Dalcin } 905*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 906*16d9e3a6SLisandro Dalcin } 907*16d9e3a6SLisandro Dalcin 908*16d9e3a6SLisandro Dalcin /*MC 909*16d9e3a6SLisandro Dalcin PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre 910*16d9e3a6SLisandro Dalcin 911*16d9e3a6SLisandro Dalcin Options Database Keys: 912*16d9e3a6SLisandro Dalcin + -pc_hypre_type - One of pilut, parasails, boomeramg, euclid 913*16d9e3a6SLisandro Dalcin - Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX 914*16d9e3a6SLisandro Dalcin preconditioner 915*16d9e3a6SLisandro Dalcin 916*16d9e3a6SLisandro Dalcin Level: intermediate 917*16d9e3a6SLisandro Dalcin 918*16d9e3a6SLisandro Dalcin Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()), 919*16d9e3a6SLisandro Dalcin the many hypre options can ONLY be set via the options database (e.g. the command line 920*16d9e3a6SLisandro Dalcin or with PetscOptionsSetValue(), there are no functions to set them) 921*16d9e3a6SLisandro Dalcin 922*16d9e3a6SLisandro Dalcin The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations 923*16d9e3a6SLisandro Dalcin (V-cycles) that boomeramg does EACH time it is called. So for example, if it is set to 2 then 924*16d9e3a6SLisandro Dalcin 2-V-cycles are being used to define the preconditioner. -ksp_max_iter and -ksp_rtol STILL determine 925*16d9e3a6SLisandro Dalcin the total number of iterations and tolerance for the Krylov solver. For example, if 926*16d9e3a6SLisandro Dalcin -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10 then AT MOST twenty V-cycles of boomeramg 927*16d9e3a6SLisandro Dalcin will be called. 928*16d9e3a6SLisandro Dalcin 929*16d9e3a6SLisandro Dalcin If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly 930*16d9e3a6SLisandro Dalcin and use -ksp_max_it to control the number of V-cycles. 931*16d9e3a6SLisandro Dalcin (see the PETSc FAQ.html at the PETSc website under the Documentation tab). 932*16d9e3a6SLisandro Dalcin 933*16d9e3a6SLisandro Dalcin 2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option 934*16d9e3a6SLisandro Dalcin -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L. 935*16d9e3a6SLisandro Dalcin 936*16d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 937*16d9e3a6SLisandro Dalcin PCHYPRESetType() 938*16d9e3a6SLisandro Dalcin 939*16d9e3a6SLisandro Dalcin M*/ 940*16d9e3a6SLisandro Dalcin 941*16d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 942*16d9e3a6SLisandro Dalcin #undef __FUNCT__ 943*16d9e3a6SLisandro Dalcin #define __FUNCT__ "PCCreate_HYPRE" 944*16d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_HYPRE(PC pc) 945*16d9e3a6SLisandro Dalcin { 946*16d9e3a6SLisandro Dalcin PC_HYPRE *jac; 947*16d9e3a6SLisandro Dalcin PetscErrorCode ierr; 948*16d9e3a6SLisandro Dalcin 949*16d9e3a6SLisandro Dalcin PetscFunctionBegin; 950*16d9e3a6SLisandro Dalcin ierr = PetscNew(PC_HYPRE,&jac);CHKERRQ(ierr); 951*16d9e3a6SLisandro Dalcin ierr = PetscLogObjectMemory(pc,sizeof(PC_HYPRE));CHKERRQ(ierr); 952*16d9e3a6SLisandro Dalcin pc->data = jac; 953*16d9e3a6SLisandro Dalcin pc->ops->destroy = PCDestroy_HYPRE; 954*16d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE; 955*16d9e3a6SLisandro Dalcin pc->ops->setup = PCSetUp_HYPRE; 956*16d9e3a6SLisandro Dalcin pc->ops->apply = PCApply_HYPRE; 957*16d9e3a6SLisandro Dalcin jac->comm_hypre = MPI_COMM_NULL; 958*16d9e3a6SLisandro Dalcin jac->hypre_type = PETSC_NULL; 959*16d9e3a6SLisandro Dalcin /* duplicate communicator for hypre */ 960*16d9e3a6SLisandro Dalcin ierr = MPI_Comm_dup(pc->comm,&(jac->comm_hypre));CHKERRQ(ierr); 961*16d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr); 962*16d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr); 963*16d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 964*16d9e3a6SLisandro Dalcin } 965*16d9e3a6SLisandro Dalcin EXTERN_C_END 966