116d9e3a6SLisandro Dalcin #define PETSCKSP_DLL 216d9e3a6SLisandro Dalcin 316d9e3a6SLisandro Dalcin /* 416d9e3a6SLisandro Dalcin Provides an interface to the LLNL package hypre 516d9e3a6SLisandro Dalcin */ 6*0f1074feSSatish Balay 7*0f1074feSSatish Balay /* Must use hypre 2.0.0 or more recent. */ 8*0f1074feSSatish Balay 916d9e3a6SLisandro Dalcin #include "private/pcimpl.h" /*I "petscpc.h" I*/ 1016d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 1116d9e3a6SLisandro Dalcin #include "HYPRE.h" 12*0f1074feSSatish Balay #include "HYPRE_parcsr_ls.h" 13*0f1074feSSatish Balay #include "_hypre_parcsr_mv.h" 14*0f1074feSSatish Balay #include "_hypre_IJ_mv.h" 1516d9e3a6SLisandro Dalcin EXTERN_C_END 1616d9e3a6SLisandro Dalcin 1716d9e3a6SLisandro Dalcin EXTERN PetscErrorCode MatHYPRE_IJMatrixCreate(Mat,HYPRE_IJMatrix*); 1816d9e3a6SLisandro Dalcin EXTERN PetscErrorCode MatHYPRE_IJMatrixCopy(Mat,HYPRE_IJMatrix); 1916d9e3a6SLisandro Dalcin EXTERN PetscErrorCode VecHYPRE_IJVectorCreate(Vec,HYPRE_IJVector*); 2016d9e3a6SLisandro Dalcin 2116d9e3a6SLisandro Dalcin /* 2216d9e3a6SLisandro Dalcin Private context (data structure) for the preconditioner. 2316d9e3a6SLisandro Dalcin */ 2416d9e3a6SLisandro Dalcin typedef struct { 2516d9e3a6SLisandro Dalcin HYPRE_Solver hsolver; 2616d9e3a6SLisandro Dalcin HYPRE_IJMatrix ij; 2716d9e3a6SLisandro Dalcin HYPRE_IJVector b,x; 2816d9e3a6SLisandro Dalcin 2916d9e3a6SLisandro Dalcin PetscErrorCode (*destroy)(HYPRE_Solver); 3016d9e3a6SLisandro Dalcin PetscErrorCode (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector); 3116d9e3a6SLisandro Dalcin PetscErrorCode (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector); 3216d9e3a6SLisandro Dalcin 3316d9e3a6SLisandro Dalcin MPI_Comm comm_hypre; 3416d9e3a6SLisandro Dalcin char *hypre_type; 3516d9e3a6SLisandro Dalcin 3616d9e3a6SLisandro Dalcin /* options for Pilut and BoomerAMG*/ 3716d9e3a6SLisandro Dalcin int maxiter; 3816d9e3a6SLisandro Dalcin double tol; 3916d9e3a6SLisandro Dalcin 4016d9e3a6SLisandro Dalcin /* options for Pilut */ 4116d9e3a6SLisandro Dalcin int factorrowsize; 4216d9e3a6SLisandro Dalcin 4316d9e3a6SLisandro Dalcin /* options for ParaSails */ 4416d9e3a6SLisandro Dalcin int nlevels; 4516d9e3a6SLisandro Dalcin double threshhold; 4616d9e3a6SLisandro Dalcin double filter; 4716d9e3a6SLisandro Dalcin int sym; 4816d9e3a6SLisandro Dalcin double loadbal; 4916d9e3a6SLisandro Dalcin int logging; 5016d9e3a6SLisandro Dalcin int ruse; 5116d9e3a6SLisandro Dalcin int symt; 5216d9e3a6SLisandro Dalcin 5316d9e3a6SLisandro Dalcin /* options for Euclid */ 5416d9e3a6SLisandro Dalcin PetscTruth bjilu; 5516d9e3a6SLisandro Dalcin int levels; 5616d9e3a6SLisandro Dalcin 5716d9e3a6SLisandro Dalcin /* options for Euclid and BoomerAMG */ 5816d9e3a6SLisandro Dalcin PetscTruth printstatistics; 5916d9e3a6SLisandro Dalcin 6016d9e3a6SLisandro Dalcin /* options for BoomerAMG */ 6116d9e3a6SLisandro Dalcin int cycletype; 6216d9e3a6SLisandro Dalcin int maxlevels; 6316d9e3a6SLisandro Dalcin double strongthreshold; 6416d9e3a6SLisandro Dalcin double maxrowsum; 65*0f1074feSSatish Balay int gridsweeps[3]; 6616d9e3a6SLisandro Dalcin int coarsentype; 6716d9e3a6SLisandro Dalcin int measuretype; 68*0f1074feSSatish Balay int relaxtype[3]; 6916d9e3a6SLisandro Dalcin double relaxweight; 7016d9e3a6SLisandro Dalcin double outerrelaxweight; 7116d9e3a6SLisandro Dalcin int relaxorder; 7216d9e3a6SLisandro Dalcin double truncfactor; 7316d9e3a6SLisandro Dalcin PetscTruth applyrichardson; 74*0f1074feSSatish Balay int pmax; 75*0f1074feSSatish Balay int interptype; 76*0f1074feSSatish Balay int agg_nl; 77*0f1074feSSatish Balay int agg_num_paths; 7816d9e3a6SLisandro Dalcin 7916d9e3a6SLisandro Dalcin } PC_HYPRE; 8016d9e3a6SLisandro Dalcin 8116d9e3a6SLisandro Dalcin 8216d9e3a6SLisandro Dalcin #undef __FUNCT__ 8316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetUp_HYPRE" 8416d9e3a6SLisandro Dalcin static PetscErrorCode PCSetUp_HYPRE(PC pc) 8516d9e3a6SLisandro Dalcin { 8616d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 8716d9e3a6SLisandro Dalcin PetscErrorCode ierr; 8816d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 8916d9e3a6SLisandro Dalcin HYPRE_ParVector bv,xv; 9016d9e3a6SLisandro Dalcin PetscInt bs; 9116d9e3a6SLisandro Dalcin int hierr; 9216d9e3a6SLisandro Dalcin 9316d9e3a6SLisandro Dalcin PetscFunctionBegin; 9416d9e3a6SLisandro Dalcin if (!jac->hypre_type) { 9516d9e3a6SLisandro Dalcin ierr = PCHYPRESetType(pc,"pilut");CHKERRQ(ierr); 9616d9e3a6SLisandro Dalcin } 9716d9e3a6SLisandro Dalcin #if 1 9816d9e3a6SLisandro Dalcin if (!pc->setupcalled) { 9916d9e3a6SLisandro Dalcin /* create the matrix and vectors the first time through */ 10016d9e3a6SLisandro Dalcin Vec x,b; 10116d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr); 10216d9e3a6SLisandro Dalcin ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr); 10316d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr); 10416d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr); 10516d9e3a6SLisandro Dalcin ierr = VecDestroy(x);CHKERRQ(ierr); 10616d9e3a6SLisandro Dalcin ierr = VecDestroy(b);CHKERRQ(ierr); 10716d9e3a6SLisandro Dalcin } else if (pc->flag != SAME_NONZERO_PATTERN) { 10816d9e3a6SLisandro Dalcin /* rebuild the matrix from scratch */ 10916d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); 11016d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr); 11116d9e3a6SLisandro Dalcin } 11216d9e3a6SLisandro Dalcin #else 11316d9e3a6SLisandro Dalcin if (!jac->ij) { /* create the matrix the first time through */ 11416d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr); 11516d9e3a6SLisandro Dalcin } 11616d9e3a6SLisandro Dalcin if (!jac->b) { /* create the vectors the first time through */ 11716d9e3a6SLisandro Dalcin Vec x,b; 11816d9e3a6SLisandro Dalcin ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr); 11916d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr); 12016d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr); 12116d9e3a6SLisandro Dalcin ierr = VecDestroy(x);CHKERRQ(ierr); 12216d9e3a6SLisandro Dalcin ierr = VecDestroy(b);CHKERRQ(ierr); 12316d9e3a6SLisandro Dalcin } 12416d9e3a6SLisandro Dalcin #endif 12516d9e3a6SLisandro Dalcin /* special case for BoomerAMG */ 12616d9e3a6SLisandro Dalcin if (jac->setup == HYPRE_BoomerAMGSetup) { 12716d9e3a6SLisandro Dalcin ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr); 12816d9e3a6SLisandro Dalcin if (bs > 1) { 12916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetNumFunctions(jac->hsolver,bs);CHKERRQ(ierr); 13016d9e3a6SLisandro Dalcin } 13116d9e3a6SLisandro Dalcin }; 13216d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr); 13316d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 13416d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&bv);CHKERRQ(ierr); 13516d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&xv);CHKERRQ(ierr); 13616d9e3a6SLisandro Dalcin hierr = (*jac->setup)(jac->hsolver,hmat,bv,xv); 13716d9e3a6SLisandro Dalcin if (hierr) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE setup, error code %d",hierr); 13816d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 13916d9e3a6SLisandro Dalcin } 14016d9e3a6SLisandro Dalcin 14116d9e3a6SLisandro Dalcin /* 14216d9e3a6SLisandro Dalcin Replaces the address where the HYPRE vector points to its data with the address of 14316d9e3a6SLisandro Dalcin PETSc's data. Saves the old address so it can be reset when we are finished with it. 14416d9e3a6SLisandro Dalcin Allows use to get the data into a HYPRE vector without the cost of memcopies 14516d9e3a6SLisandro Dalcin */ 14616d9e3a6SLisandro Dalcin #define HYPREReplacePointer(b,newvalue,savedvalue) {\ 14716d9e3a6SLisandro Dalcin hypre_ParVector *par_vector = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));\ 14816d9e3a6SLisandro Dalcin hypre_Vector *local_vector = hypre_ParVectorLocalVector(par_vector);\ 14916d9e3a6SLisandro Dalcin savedvalue = local_vector->data;\ 15016d9e3a6SLisandro Dalcin local_vector->data = newvalue;} 15116d9e3a6SLisandro Dalcin 15216d9e3a6SLisandro Dalcin #undef __FUNCT__ 15316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApply_HYPRE" 15416d9e3a6SLisandro Dalcin static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x) 15516d9e3a6SLisandro Dalcin { 15616d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 15716d9e3a6SLisandro Dalcin PetscErrorCode ierr; 15816d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 15916d9e3a6SLisandro Dalcin PetscScalar *bv,*xv; 16016d9e3a6SLisandro Dalcin HYPRE_ParVector jbv,jxv; 16116d9e3a6SLisandro Dalcin PetscScalar *sbv,*sxv; 16216d9e3a6SLisandro Dalcin int hierr; 16316d9e3a6SLisandro Dalcin 16416d9e3a6SLisandro Dalcin PetscFunctionBegin; 16516d9e3a6SLisandro Dalcin if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);} 16616d9e3a6SLisandro Dalcin ierr = VecGetArray(b,&bv);CHKERRQ(ierr); 16716d9e3a6SLisandro Dalcin ierr = VecGetArray(x,&xv);CHKERRQ(ierr); 16816d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,bv,sbv); 16916d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,xv,sxv); 17016d9e3a6SLisandro Dalcin 17116d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 17216d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr); 17316d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr); 17416d9e3a6SLisandro Dalcin hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv); 175*0f1074feSSatish Balay 176*0f1074feSSatish Balay /*if (hierr && (hierr != HYPRE_ERROR_CONV || jac->solve != HYPRE_BoomerAMGSolve))SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 177*0f1074feSSatish Balay */ 178*0f1074feSSatish Balay /* error code of HYPRE_ERROR_CONV means convergence not achieved - if 179*0f1074feSSatish Balay the tolerance is set to 0.0 (the default), a convergence error will 180*0f1074feSSatish Balay not occur (so we may not want to overide the conv. error here?*/ 181*0f1074feSSatish Balay if (hierr && hierr != HYPRE_ERROR_CONV) 182*0f1074feSSatish Balay { 183*0f1074feSSatish Balay SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 184*0f1074feSSatish Balay } 18516d9e3a6SLisandro Dalcin if (hierr) hypre__global_error = 0; 18616d9e3a6SLisandro Dalcin 187*0f1074feSSatish Balay 18816d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,sbv,bv); 18916d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,sxv,xv); 19016d9e3a6SLisandro Dalcin ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr); 19116d9e3a6SLisandro Dalcin ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr); 19216d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 19316d9e3a6SLisandro Dalcin } 19416d9e3a6SLisandro Dalcin 19516d9e3a6SLisandro Dalcin #undef __FUNCT__ 19616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCDestroy_HYPRE" 19716d9e3a6SLisandro Dalcin static PetscErrorCode PCDestroy_HYPRE(PC pc) 19816d9e3a6SLisandro Dalcin { 19916d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 20016d9e3a6SLisandro Dalcin PetscErrorCode ierr; 20116d9e3a6SLisandro Dalcin 20216d9e3a6SLisandro Dalcin PetscFunctionBegin; 20316d9e3a6SLisandro Dalcin if (jac->ij) { ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); } 20416d9e3a6SLisandro Dalcin if (jac->b) { ierr = HYPRE_IJVectorDestroy(jac->b);CHKERRQ(ierr); } 20516d9e3a6SLisandro Dalcin if (jac->x) { ierr = HYPRE_IJVectorDestroy(jac->x);CHKERRQ(ierr); } 20616d9e3a6SLisandro Dalcin if (jac->destroy) { ierr = (*jac->destroy)(jac->hsolver);CHKERRQ(ierr); } 20716d9e3a6SLisandro Dalcin ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr); 20816d9e3a6SLisandro Dalcin if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);} 20916d9e3a6SLisandro Dalcin ierr = PetscFree(jac);CHKERRQ(ierr); 21016d9e3a6SLisandro Dalcin 21116d9e3a6SLisandro Dalcin ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr); 21216d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);CHKERRQ(ierr); 21316d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);CHKERRQ(ierr); 21416d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 21516d9e3a6SLisandro Dalcin } 21616d9e3a6SLisandro Dalcin 21716d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 21816d9e3a6SLisandro Dalcin #undef __FUNCT__ 21916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut" 22016d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc) 22116d9e3a6SLisandro Dalcin { 22216d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 22316d9e3a6SLisandro Dalcin PetscErrorCode ierr; 22416d9e3a6SLisandro Dalcin PetscTruth flag; 22516d9e3a6SLisandro Dalcin 22616d9e3a6SLisandro Dalcin PetscFunctionBegin; 22716d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr); 22816d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr); 22916d9e3a6SLisandro Dalcin if (flag) { 23016d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 23116d9e3a6SLisandro Dalcin } 23216d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr); 23316d9e3a6SLisandro Dalcin if (flag) { 23416d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutSetDropTolerance(jac->hsolver,jac->tol);CHKERRQ(ierr); 23516d9e3a6SLisandro Dalcin } 23616d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr); 23716d9e3a6SLisandro Dalcin if (flag) { 23816d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutSetFactorRowSize(jac->hsolver,jac->factorrowsize);CHKERRQ(ierr); 23916d9e3a6SLisandro Dalcin } 24016d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 24116d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 24216d9e3a6SLisandro Dalcin } 24316d9e3a6SLisandro Dalcin 24416d9e3a6SLisandro Dalcin #undef __FUNCT__ 24516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Pilut" 24616d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer) 24716d9e3a6SLisandro Dalcin { 24816d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 24916d9e3a6SLisandro Dalcin PetscErrorCode ierr; 25016d9e3a6SLisandro Dalcin PetscTruth iascii; 25116d9e3a6SLisandro Dalcin 25216d9e3a6SLisandro Dalcin PetscFunctionBegin; 25316d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 25416d9e3a6SLisandro Dalcin if (iascii) { 25516d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut preconditioning\n");CHKERRQ(ierr); 25616d9e3a6SLisandro Dalcin if (jac->maxiter != PETSC_DEFAULT) { 25716d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr); 25816d9e3a6SLisandro Dalcin } else { 25916d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr); 26016d9e3a6SLisandro Dalcin } 26116d9e3a6SLisandro Dalcin if (jac->tol != PETSC_DEFAULT) { 26216d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr); 26316d9e3a6SLisandro Dalcin } else { 26416d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr); 26516d9e3a6SLisandro Dalcin } 26616d9e3a6SLisandro Dalcin if (jac->factorrowsize != PETSC_DEFAULT) { 26716d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr); 26816d9e3a6SLisandro Dalcin } else { 26916d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default factor row size \n");CHKERRQ(ierr); 27016d9e3a6SLisandro Dalcin } 27116d9e3a6SLisandro Dalcin } 27216d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 27316d9e3a6SLisandro Dalcin } 27416d9e3a6SLisandro Dalcin 27516d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 27616d9e3a6SLisandro Dalcin #undef __FUNCT__ 27716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid" 27816d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc) 27916d9e3a6SLisandro Dalcin { 28016d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 28116d9e3a6SLisandro Dalcin PetscErrorCode ierr; 28216d9e3a6SLisandro Dalcin PetscTruth flag; 283390e7148SBarry Smith char *args[8],levels[16]; 284390e7148SBarry Smith PetscInt cnt = 0; 28516d9e3a6SLisandro Dalcin 28616d9e3a6SLisandro Dalcin PetscFunctionBegin; 28716d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr); 28816d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr); 28916d9e3a6SLisandro Dalcin if (flag) { 29016d9e3a6SLisandro Dalcin if (jac->levels < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels); 29116d9e3a6SLisandro Dalcin sprintf(levels,"%d",jac->levels); 292390e7148SBarry Smith args[cnt++] = (char*)"-level"; args[cnt++] = levels; 29316d9e3a6SLisandro Dalcin } 29416d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);CHKERRQ(ierr); 29516d9e3a6SLisandro Dalcin if (jac->bjilu) { 296390e7148SBarry Smith args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1"; 29716d9e3a6SLisandro Dalcin } 29816d9e3a6SLisandro Dalcin 29916d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);CHKERRQ(ierr); 30016d9e3a6SLisandro Dalcin if (jac->printstatistics) { 301390e7148SBarry Smith args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1"; 302390e7148SBarry Smith args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1"; 30316d9e3a6SLisandro Dalcin } 30416d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 305390e7148SBarry Smith if (cnt) { 306390e7148SBarry Smith ierr = HYPRE_EuclidSetParams(jac->hsolver,cnt,args);CHKERRQ(ierr); 307390e7148SBarry Smith } 30816d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 30916d9e3a6SLisandro Dalcin } 31016d9e3a6SLisandro Dalcin 31116d9e3a6SLisandro Dalcin #undef __FUNCT__ 31216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Euclid" 31316d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer) 31416d9e3a6SLisandro Dalcin { 31516d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 31616d9e3a6SLisandro Dalcin PetscErrorCode ierr; 31716d9e3a6SLisandro Dalcin PetscTruth iascii; 31816d9e3a6SLisandro Dalcin 31916d9e3a6SLisandro Dalcin PetscFunctionBegin; 32016d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 32116d9e3a6SLisandro Dalcin if (iascii) { 32216d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid preconditioning\n");CHKERRQ(ierr); 32316d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr); 32416d9e3a6SLisandro Dalcin if (jac->bjilu) { 32516d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr); 32616d9e3a6SLisandro Dalcin } 32716d9e3a6SLisandro Dalcin } 32816d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 32916d9e3a6SLisandro Dalcin } 33016d9e3a6SLisandro Dalcin 33116d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 33216d9e3a6SLisandro Dalcin 33316d9e3a6SLisandro Dalcin #undef __FUNCT__ 33416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG" 33516d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x) 33616d9e3a6SLisandro Dalcin { 33716d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 33816d9e3a6SLisandro Dalcin PetscErrorCode ierr; 33916d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 34016d9e3a6SLisandro Dalcin PetscScalar *bv,*xv; 34116d9e3a6SLisandro Dalcin HYPRE_ParVector jbv,jxv; 34216d9e3a6SLisandro Dalcin PetscScalar *sbv,*sxv; 34316d9e3a6SLisandro Dalcin int hierr; 34416d9e3a6SLisandro Dalcin 34516d9e3a6SLisandro Dalcin PetscFunctionBegin; 34616d9e3a6SLisandro Dalcin ierr = VecSet(x,0.0);CHKERRQ(ierr); 34716d9e3a6SLisandro Dalcin ierr = VecGetArray(b,&bv);CHKERRQ(ierr); 34816d9e3a6SLisandro Dalcin ierr = VecGetArray(x,&xv);CHKERRQ(ierr); 34916d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,bv,sbv); 35016d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,xv,sxv); 35116d9e3a6SLisandro Dalcin 35216d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 35316d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr); 35416d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr); 35516d9e3a6SLisandro Dalcin 35616d9e3a6SLisandro Dalcin hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv); 35716d9e3a6SLisandro Dalcin /* error code of 1 in BoomerAMG merely means convergence not achieved */ 35816d9e3a6SLisandro Dalcin if (hierr && (hierr != 1)) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 35916d9e3a6SLisandro Dalcin if (hierr) hypre__global_error = 0; 36016d9e3a6SLisandro Dalcin 36116d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,sbv,bv); 36216d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,sxv,xv); 36316d9e3a6SLisandro Dalcin ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr); 36416d9e3a6SLisandro Dalcin ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr); 36516d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 36616d9e3a6SLisandro Dalcin } 36716d9e3a6SLisandro Dalcin 36816d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCycleType[] = {"","V","W"}; 369*0f1074feSSatish Balay static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"}; 37016d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGMeasureType[] = {"local","global"}; 371*0f1074feSSatish Balay static const char *HYPREBoomerAMGRelaxType[] = {"Jacobi","sequential-Gauss-Seidel","","SOR/Jacobi","backward-SOR/Jacobi","","symmetric-SOR/Jacobi", 372*0f1074feSSatish Balay "","","Gaussian-elimination"}; 373*0f1074feSSatish Balay static const char *HYPREBoomerAMGInterpType[] = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i", 374*0f1074feSSatish Balay "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"}; 37516d9e3a6SLisandro Dalcin #undef __FUNCT__ 37616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG" 37716d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc) 37816d9e3a6SLisandro Dalcin { 37916d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 38016d9e3a6SLisandro Dalcin PetscErrorCode ierr; 38116d9e3a6SLisandro Dalcin int n,indx; 38216d9e3a6SLisandro Dalcin PetscTruth flg, tmp_truth; 38316d9e3a6SLisandro Dalcin double tmpdbl, twodbl[2]; 38416d9e3a6SLisandro Dalcin 38516d9e3a6SLisandro Dalcin PetscFunctionBegin; 38616d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr); 38716d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr); 38816d9e3a6SLisandro Dalcin if (flg) { 38916d9e3a6SLisandro Dalcin jac->cycletype = indx; 39016d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr); 39116d9e3a6SLisandro Dalcin } 39216d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr); 39316d9e3a6SLisandro Dalcin if (flg) { 39416d9e3a6SLisandro Dalcin if (jac->maxlevels < 2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels); 39516d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr); 39616d9e3a6SLisandro Dalcin } 39716d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr); 39816d9e3a6SLisandro Dalcin if (flg) { 39916d9e3a6SLisandro Dalcin if (jac->maxiter < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter); 40016d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 40116d9e3a6SLisandro Dalcin } 402*0f1074feSSatish 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); 40316d9e3a6SLisandro Dalcin if (flg) { 404*0f1074feSSatish Balay if (jac->tol < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be greater than or equal to zero",jac->tol); 40516d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 40616d9e3a6SLisandro Dalcin } 40716d9e3a6SLisandro Dalcin 408*0f1074feSSatish Balay ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr); 40916d9e3a6SLisandro Dalcin if (flg) { 41016d9e3a6SLisandro Dalcin if (jac->truncfactor < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor); 41116d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr); 41216d9e3a6SLisandro Dalcin } 41316d9e3a6SLisandro Dalcin 414*0f1074feSSatish 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); 415*0f1074feSSatish Balay if (flg) { 416*0f1074feSSatish Balay if (jac->pmax < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"P_max %G must be greater than or equal to zero",jac->pmax); 417*0f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);CHKERRQ(ierr); 418*0f1074feSSatish Balay } 419*0f1074feSSatish Balay 420*0f1074feSSatish Balay ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);CHKERRQ(ierr); 421*0f1074feSSatish Balay if (flg) { 422*0f1074feSSatish Balay if (jac->agg_nl < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %G must be greater than or equal to zero",jac->agg_nl); 423*0f1074feSSatish Balay 424*0f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl);CHKERRQ(ierr); 425*0f1074feSSatish Balay } 426*0f1074feSSatish Balay 427*0f1074feSSatish Balay 428*0f1074feSSatish 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); 429*0f1074feSSatish Balay if (flg) { 430*0f1074feSSatish Balay if (jac->agg_num_paths < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of paths %G must be greater than or equal to 1",jac->agg_num_paths); 431*0f1074feSSatish Balay 432*0f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);CHKERRQ(ierr); 433*0f1074feSSatish Balay } 434*0f1074feSSatish Balay 435*0f1074feSSatish Balay 43616d9e3a6SLisandro Dalcin ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr); 43716d9e3a6SLisandro Dalcin if (flg) { 43816d9e3a6SLisandro Dalcin if (jac->strongthreshold < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold); 43916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr); 44016d9e3a6SLisandro Dalcin } 44116d9e3a6SLisandro Dalcin ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr); 44216d9e3a6SLisandro Dalcin if (flg) { 44316d9e3a6SLisandro Dalcin if (jac->maxrowsum < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum); 44416d9e3a6SLisandro Dalcin if (jac->maxrowsum > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum); 44516d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr); 44616d9e3a6SLisandro Dalcin } 44716d9e3a6SLisandro Dalcin 44816d9e3a6SLisandro Dalcin /* Grid sweeps */ 449*0f1074feSSatish 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); 45016d9e3a6SLisandro Dalcin if (flg) { 45116d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver,indx);CHKERRQ(ierr); 45216d9e3a6SLisandro Dalcin /* modify the jac structure so we can view the updated options with PC_View */ 45316d9e3a6SLisandro Dalcin jac->gridsweeps[0] = indx; 454*0f1074feSSatish Balay jac->gridsweeps[1] = indx; 455*0f1074feSSatish Balay /*defaults coarse to 1 */ 456*0f1074feSSatish Balay jac->gridsweeps[2] = 1; 45716d9e3a6SLisandro Dalcin } 458*0f1074feSSatish Balay 459*0f1074feSSatish Balay ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx ,&flg);CHKERRQ(ierr); 46016d9e3a6SLisandro Dalcin if (flg) { 46116d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 1);CHKERRQ(ierr); 462*0f1074feSSatish Balay jac->gridsweeps[0] = indx; 46316d9e3a6SLisandro Dalcin } 46416d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr); 46516d9e3a6SLisandro Dalcin if (flg) { 46616d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 2);CHKERRQ(ierr); 467*0f1074feSSatish Balay jac->gridsweeps[1] = indx; 46816d9e3a6SLisandro Dalcin } 469*0f1074feSSatish Balay ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr); 47016d9e3a6SLisandro Dalcin if (flg) { 47116d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 3);CHKERRQ(ierr); 472*0f1074feSSatish Balay jac->gridsweeps[2] = indx; 47316d9e3a6SLisandro Dalcin } 47416d9e3a6SLisandro Dalcin 47516d9e3a6SLisandro Dalcin /* Relax type */ 476*0f1074feSSatish 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); 47716d9e3a6SLisandro Dalcin if (flg) { 478*0f1074feSSatish Balay jac->relaxtype[0] = jac->relaxtype[1] = indx; 47916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, indx);CHKERRQ(ierr); 480*0f1074feSSatish Balay /* by default, coarse type set to 9 */ 481*0f1074feSSatish Balay jac->relaxtype[2] = 9; 482*0f1074feSSatish Balay 48316d9e3a6SLisandro Dalcin } 484*0f1074feSSatish Balay ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr); 48516d9e3a6SLisandro Dalcin if (flg) { 48616d9e3a6SLisandro Dalcin jac->relaxtype[0] = indx; 48716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 1);CHKERRQ(ierr); 48816d9e3a6SLisandro Dalcin } 489*0f1074feSSatish Balay ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr); 49016d9e3a6SLisandro Dalcin if (flg) { 491*0f1074feSSatish Balay jac->relaxtype[1] = indx; 49216d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 2);CHKERRQ(ierr); 49316d9e3a6SLisandro Dalcin } 49416d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr); 49516d9e3a6SLisandro Dalcin if (flg) { 496*0f1074feSSatish Balay jac->relaxtype[2] = indx; 49716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 3);CHKERRQ(ierr); 49816d9e3a6SLisandro Dalcin } 49916d9e3a6SLisandro Dalcin 50016d9e3a6SLisandro Dalcin /* Relaxation Weight */ 50116d9e3a6SLisandro 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); 50216d9e3a6SLisandro Dalcin if (flg) { 50316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxWt(jac->hsolver,tmpdbl);CHKERRQ(ierr); 50416d9e3a6SLisandro Dalcin jac->relaxweight = tmpdbl; 50516d9e3a6SLisandro Dalcin } 50616d9e3a6SLisandro Dalcin 50716d9e3a6SLisandro Dalcin n=2; 50816d9e3a6SLisandro Dalcin twodbl[0] = twodbl[1] = 1.0; 50916d9e3a6SLisandro 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); 51016d9e3a6SLisandro Dalcin if (flg) { 51116d9e3a6SLisandro Dalcin if (n == 2) { 51216d9e3a6SLisandro Dalcin indx = (int)PetscAbsReal(twodbl[1]); 51316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetLevelRelaxWt(jac->hsolver,twodbl[0],indx);CHKERRQ(ierr); 51416d9e3a6SLisandro Dalcin } else { 51516d9e3a6SLisandro 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); 51616d9e3a6SLisandro Dalcin } 51716d9e3a6SLisandro Dalcin } 51816d9e3a6SLisandro Dalcin 51916d9e3a6SLisandro Dalcin /* Outer relaxation Weight */ 52016d9e3a6SLisandro 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); 52116d9e3a6SLisandro Dalcin if (flg) { 52216d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetOuterWt( jac->hsolver, tmpdbl);CHKERRQ(ierr); 52316d9e3a6SLisandro Dalcin jac->outerrelaxweight = tmpdbl; 52416d9e3a6SLisandro Dalcin } 52516d9e3a6SLisandro Dalcin 52616d9e3a6SLisandro Dalcin n=2; 52716d9e3a6SLisandro Dalcin twodbl[0] = twodbl[1] = 1.0; 52816d9e3a6SLisandro 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); 52916d9e3a6SLisandro Dalcin if (flg) { 53016d9e3a6SLisandro Dalcin if (n == 2) { 53116d9e3a6SLisandro Dalcin indx = (int)PetscAbsReal(twodbl[1]); 53216d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetLevelOuterWt( jac->hsolver, twodbl[0], indx);CHKERRQ(ierr); 53316d9e3a6SLisandro Dalcin } else { 53416d9e3a6SLisandro 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); 53516d9e3a6SLisandro Dalcin } 53616d9e3a6SLisandro Dalcin } 53716d9e3a6SLisandro Dalcin 53816d9e3a6SLisandro Dalcin /* the Relax Order */ 53916d9e3a6SLisandro Dalcin /* ierr = PetscOptionsName("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", &flg);CHKERRQ(ierr); */ 54016d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr); 54116d9e3a6SLisandro Dalcin 54216d9e3a6SLisandro Dalcin if (flg) { 54316d9e3a6SLisandro Dalcin jac->relaxorder = 0; 54416d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr); 54516d9e3a6SLisandro Dalcin } 54616d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr); 54716d9e3a6SLisandro Dalcin if (flg) { 54816d9e3a6SLisandro Dalcin jac->measuretype = indx; 54916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr); 55016d9e3a6SLisandro Dalcin } 551*0f1074feSSatish Balay /* update list length 3/07 */ 552*0f1074feSSatish Balay ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,11,HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr); 55316d9e3a6SLisandro Dalcin if (flg) { 55416d9e3a6SLisandro Dalcin jac->coarsentype = indx; 55516d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr); 55616d9e3a6SLisandro Dalcin } 557*0f1074feSSatish Balay 558*0f1074feSSatish Balay /* new 3/07 */ 559*0f1074feSSatish Balay ierr = PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,14,HYPREBoomerAMGInterpType[0],&indx,&flg);CHKERRQ(ierr); 560*0f1074feSSatish Balay if (flg) { 561*0f1074feSSatish Balay jac->interptype = indx; 562*0f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);CHKERRQ(ierr); 563*0f1074feSSatish Balay } 564*0f1074feSSatish Balay 565*0f1074feSSatish Balay 566*0f1074feSSatish Balay 56716d9e3a6SLisandro Dalcin ierr = PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);CHKERRQ(ierr); 56816d9e3a6SLisandro Dalcin if (flg) { 56916d9e3a6SLisandro Dalcin int level=3; 57016d9e3a6SLisandro Dalcin jac->printstatistics = PETSC_TRUE; 57116d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);CHKERRQ(ierr); 57216d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetPrintLevel(jac->hsolver,level);CHKERRQ(ierr); 57316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetDebugFlag(jac->hsolver,level);CHKERRQ(ierr); 57416d9e3a6SLisandro Dalcin } 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" 58116d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its) 58216d9e3a6SLisandro Dalcin { 58316d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 58416d9e3a6SLisandro Dalcin PetscErrorCode ierr; 58516d9e3a6SLisandro Dalcin 58616d9e3a6SLisandro Dalcin PetscFunctionBegin; 58716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,its*jac->maxiter);CHKERRQ(ierr); 58816d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,PetscMin(rtol,jac->tol));CHKERRQ(ierr); 58916d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_TRUE; 59016d9e3a6SLisandro Dalcin ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr); 59116d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_FALSE; 59216d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 59316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 59416d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 59516d9e3a6SLisandro Dalcin } 59616d9e3a6SLisandro Dalcin 59716d9e3a6SLisandro Dalcin 59816d9e3a6SLisandro Dalcin #undef __FUNCT__ 59916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_BoomerAMG" 60016d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer) 60116d9e3a6SLisandro Dalcin { 60216d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 60316d9e3a6SLisandro Dalcin PetscErrorCode ierr; 60416d9e3a6SLisandro Dalcin PetscTruth iascii; 60516d9e3a6SLisandro Dalcin 60616d9e3a6SLisandro Dalcin PetscFunctionBegin; 60716d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 60816d9e3a6SLisandro Dalcin if (iascii) { 60916d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr); 61016d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr); 61116d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr); 61216d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr); 61316d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr); 61416d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr); 615*0f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation truncation factor %G\n",jac->truncfactor);CHKERRQ(ierr); 616*0f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);CHKERRQ(ierr); 617*0f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);CHKERRQ(ierr); 618*0f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);CHKERRQ(ierr); 619*0f1074feSSatish Balay 62016d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr); 62116d9e3a6SLisandro Dalcin 622*0f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps down %d\n",jac->gridsweeps[0]);CHKERRQ(ierr); 623*0f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps up %d\n",jac->gridsweeps[1]);CHKERRQ(ierr); 624*0f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps on coarse %d\n",jac->gridsweeps[2]);CHKERRQ(ierr); 62516d9e3a6SLisandro Dalcin 626*0f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax down %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr); 627*0f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax up %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr); 628*0f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax on coarse %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr); 62916d9e3a6SLisandro Dalcin 63016d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax weight (all) %G\n",jac->relaxweight);CHKERRQ(ierr); 63116d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr); 63216d9e3a6SLisandro Dalcin 63316d9e3a6SLisandro Dalcin if (jac->relaxorder) { 63416d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr); 63516d9e3a6SLisandro Dalcin } else { 63616d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr); 63716d9e3a6SLisandro Dalcin } 63816d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Measure type %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr); 63916d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Coarsen type %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr); 640*0f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation type %s\n",HYPREBoomerAMGInterpType[jac->interptype]);CHKERRQ(ierr); 641*0f1074feSSatish Balay 64216d9e3a6SLisandro Dalcin } 64316d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 64416d9e3a6SLisandro Dalcin } 64516d9e3a6SLisandro Dalcin 64616d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 64716d9e3a6SLisandro Dalcin #undef __FUNCT__ 64816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails" 64916d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc) 65016d9e3a6SLisandro Dalcin { 65116d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 65216d9e3a6SLisandro Dalcin PetscErrorCode ierr; 65316d9e3a6SLisandro Dalcin int indx; 65416d9e3a6SLisandro Dalcin PetscTruth flag; 65516d9e3a6SLisandro Dalcin const char *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"}; 65616d9e3a6SLisandro Dalcin 65716d9e3a6SLisandro Dalcin PetscFunctionBegin; 65816d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr); 65916d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr); 66016d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr); 66116d9e3a6SLisandro Dalcin if (flag) { 66216d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr); 66316d9e3a6SLisandro Dalcin } 66416d9e3a6SLisandro Dalcin 66516d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr); 66616d9e3a6SLisandro Dalcin if (flag) { 66716d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr); 66816d9e3a6SLisandro Dalcin } 66916d9e3a6SLisandro Dalcin 67016d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr); 67116d9e3a6SLisandro Dalcin if (flag) { 67216d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr); 67316d9e3a6SLisandro Dalcin } 67416d9e3a6SLisandro Dalcin 67516d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,&flag);CHKERRQ(ierr); 67616d9e3a6SLisandro Dalcin if (flag) { 67716d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr); 67816d9e3a6SLisandro Dalcin } 67916d9e3a6SLisandro Dalcin 68016d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,&flag);CHKERRQ(ierr); 68116d9e3a6SLisandro Dalcin if (flag) { 68216d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr); 68316d9e3a6SLisandro Dalcin } 68416d9e3a6SLisandro Dalcin 68516d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);CHKERRQ(ierr); 68616d9e3a6SLisandro Dalcin if (flag) { 68716d9e3a6SLisandro Dalcin jac->symt = indx; 68816d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr); 68916d9e3a6SLisandro Dalcin } 69016d9e3a6SLisandro Dalcin 69116d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 69216d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 69316d9e3a6SLisandro Dalcin } 69416d9e3a6SLisandro Dalcin 69516d9e3a6SLisandro Dalcin #undef __FUNCT__ 69616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_ParaSails" 69716d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer) 69816d9e3a6SLisandro Dalcin { 69916d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 70016d9e3a6SLisandro Dalcin PetscErrorCode ierr; 70116d9e3a6SLisandro Dalcin PetscTruth iascii; 70216d9e3a6SLisandro Dalcin const char *symt = 0;; 70316d9e3a6SLisandro Dalcin 70416d9e3a6SLisandro Dalcin PetscFunctionBegin; 70516d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 70616d9e3a6SLisandro Dalcin if (iascii) { 70716d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails preconditioning\n");CHKERRQ(ierr); 70816d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr); 70916d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr); 71016d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr); 71116d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr); 71216d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: reuse nonzero structure %s\n",PetscTruths[jac->ruse]);CHKERRQ(ierr); 71316d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: print info to screen %s\n",PetscTruths[jac->logging]);CHKERRQ(ierr); 71416d9e3a6SLisandro Dalcin if (!jac->symt) { 71516d9e3a6SLisandro Dalcin symt = "nonsymmetric matrix and preconditioner"; 71616d9e3a6SLisandro Dalcin } else if (jac->symt == 1) { 71716d9e3a6SLisandro Dalcin symt = "SPD matrix and preconditioner"; 71816d9e3a6SLisandro Dalcin } else if (jac->symt == 2) { 71916d9e3a6SLisandro Dalcin symt = "nonsymmetric matrix but SPD preconditioner"; 72016d9e3a6SLisandro Dalcin } else { 72116d9e3a6SLisandro Dalcin SETERRQ1(PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt); 72216d9e3a6SLisandro Dalcin } 72316d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr); 72416d9e3a6SLisandro Dalcin } 72516d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 72616d9e3a6SLisandro Dalcin } 72716d9e3a6SLisandro Dalcin /* ---------------------------------------------------------------------------------*/ 72816d9e3a6SLisandro Dalcin 72916d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 73016d9e3a6SLisandro Dalcin #undef __FUNCT__ 73116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType_HYPRE" 73216d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType_HYPRE(PC pc,const char *name[]) 73316d9e3a6SLisandro Dalcin { 73416d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 73516d9e3a6SLisandro Dalcin 73616d9e3a6SLisandro Dalcin PetscFunctionBegin; 73716d9e3a6SLisandro Dalcin *name = jac->hypre_type; 73816d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 73916d9e3a6SLisandro Dalcin } 74016d9e3a6SLisandro Dalcin EXTERN_C_END 74116d9e3a6SLisandro Dalcin 74216d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 74316d9e3a6SLisandro Dalcin #undef __FUNCT__ 74416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType_HYPRE" 74516d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType_HYPRE(PC pc,const char name[]) 74616d9e3a6SLisandro Dalcin { 74716d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 74816d9e3a6SLisandro Dalcin PetscErrorCode ierr; 74916d9e3a6SLisandro Dalcin PetscTruth flag; 75016d9e3a6SLisandro Dalcin 75116d9e3a6SLisandro Dalcin PetscFunctionBegin; 75216d9e3a6SLisandro Dalcin if (jac->hypre_type) { 75316d9e3a6SLisandro Dalcin ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr); 75416d9e3a6SLisandro Dalcin if (!flag) SETERRQ(PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set"); 75516d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 75616d9e3a6SLisandro Dalcin } else { 75716d9e3a6SLisandro Dalcin ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr); 75816d9e3a6SLisandro Dalcin } 75916d9e3a6SLisandro Dalcin 76016d9e3a6SLisandro Dalcin jac->maxiter = PETSC_DEFAULT; 76116d9e3a6SLisandro Dalcin jac->tol = PETSC_DEFAULT; 76216d9e3a6SLisandro Dalcin jac->printstatistics = PetscLogPrintInfo; 76316d9e3a6SLisandro Dalcin 76416d9e3a6SLisandro Dalcin ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr); 76516d9e3a6SLisandro Dalcin if (flag) { 76616d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutCreate(jac->comm_hypre,&jac->hsolver); 76716d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut; 76816d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_Pilut; 76916d9e3a6SLisandro Dalcin jac->destroy = HYPRE_ParCSRPilutDestroy; 77016d9e3a6SLisandro Dalcin jac->setup = HYPRE_ParCSRPilutSetup; 77116d9e3a6SLisandro Dalcin jac->solve = HYPRE_ParCSRPilutSolve; 77216d9e3a6SLisandro Dalcin jac->factorrowsize = PETSC_DEFAULT; 77316d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 77416d9e3a6SLisandro Dalcin } 77516d9e3a6SLisandro Dalcin ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr); 77616d9e3a6SLisandro Dalcin if (flag) { 77716d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsCreate(jac->comm_hypre,&jac->hsolver); 77816d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails; 77916d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_ParaSails; 78016d9e3a6SLisandro Dalcin jac->destroy = HYPRE_ParaSailsDestroy; 78116d9e3a6SLisandro Dalcin jac->setup = HYPRE_ParaSailsSetup; 78216d9e3a6SLisandro Dalcin jac->solve = HYPRE_ParaSailsSolve; 78316d9e3a6SLisandro Dalcin /* initialize */ 78416d9e3a6SLisandro Dalcin jac->nlevels = 1; 78516d9e3a6SLisandro Dalcin jac->threshhold = .1; 78616d9e3a6SLisandro Dalcin jac->filter = .1; 78716d9e3a6SLisandro Dalcin jac->loadbal = 0; 78816d9e3a6SLisandro Dalcin if (PetscLogPrintInfo) { 78916d9e3a6SLisandro Dalcin jac->logging = (int) PETSC_TRUE; 79016d9e3a6SLisandro Dalcin } else { 79116d9e3a6SLisandro Dalcin jac->logging = (int) PETSC_FALSE; 79216d9e3a6SLisandro Dalcin } 79316d9e3a6SLisandro Dalcin jac->ruse = (int) PETSC_FALSE; 79416d9e3a6SLisandro Dalcin jac->symt = 0; 79516d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr); 79616d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr); 79716d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr); 79816d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr); 79916d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr); 80016d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr); 80116d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 80216d9e3a6SLisandro Dalcin } 80316d9e3a6SLisandro Dalcin ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr); 80416d9e3a6SLisandro Dalcin if (flag) { 80516d9e3a6SLisandro Dalcin ierr = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver); 80616d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid; 80716d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_Euclid; 80816d9e3a6SLisandro Dalcin jac->destroy = HYPRE_EuclidDestroy; 80916d9e3a6SLisandro Dalcin jac->setup = HYPRE_EuclidSetup; 81016d9e3a6SLisandro Dalcin jac->solve = HYPRE_EuclidSolve; 81116d9e3a6SLisandro Dalcin /* initialization */ 81216d9e3a6SLisandro Dalcin jac->bjilu = PETSC_FALSE; 81316d9e3a6SLisandro Dalcin jac->levels = 1; 81416d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 81516d9e3a6SLisandro Dalcin } 81616d9e3a6SLisandro Dalcin ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr); 81716d9e3a6SLisandro Dalcin if (flag) { 81816d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGCreate(&jac->hsolver); 81916d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_BoomerAMG; 82016d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_BoomerAMG; 82116d9e3a6SLisandro Dalcin pc->ops->applytranspose = PCApplyTranspose_HYPRE_BoomerAMG; 82216d9e3a6SLisandro Dalcin pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG; 82316d9e3a6SLisandro Dalcin jac->destroy = HYPRE_BoomerAMGDestroy; 82416d9e3a6SLisandro Dalcin jac->setup = HYPRE_BoomerAMGSetup; 82516d9e3a6SLisandro Dalcin jac->solve = HYPRE_BoomerAMGSolve; 82616d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_FALSE; 82716d9e3a6SLisandro Dalcin /* these defaults match the hypre defaults */ 82816d9e3a6SLisandro Dalcin jac->cycletype = 1; 82916d9e3a6SLisandro Dalcin jac->maxlevels = 25; 83016d9e3a6SLisandro Dalcin jac->maxiter = 1; 831*0f1074feSSatish Balay jac->tol = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses 832*0f1074feSSatish Balay convergence errors) */ 83316d9e3a6SLisandro Dalcin jac->truncfactor = 0.0; 83416d9e3a6SLisandro Dalcin jac->strongthreshold = .25; 83516d9e3a6SLisandro Dalcin jac->maxrowsum = .9; 83616d9e3a6SLisandro Dalcin jac->coarsentype = 6; 83716d9e3a6SLisandro Dalcin jac->measuretype = 0; 838*0f1074feSSatish Balay jac->gridsweeps[0] = jac->gridsweeps[1] = jac->gridsweeps[2] = 1; 839*0f1074feSSatish Balay jac->relaxtype[0] = jac->relaxtype[1] = 6; /* Now defaults to SYMMETRIC since 840*0f1074feSSatish Balay in PETSc we are using a a PC - most likely 841*0f1074feSSatish Balay with CG */ 842*0f1074feSSatish Balay jac->relaxtype[2] = 9; /*G.E. */ 84316d9e3a6SLisandro Dalcin jac->relaxweight = 1.0; 84416d9e3a6SLisandro Dalcin jac->outerrelaxweight = 1.0; 84516d9e3a6SLisandro Dalcin jac->relaxorder = 1; 846*0f1074feSSatish Balay jac->interptype = 0; 847*0f1074feSSatish Balay jac->agg_nl = 0; 848*0f1074feSSatish Balay jac->pmax = 0; 849*0f1074feSSatish Balay jac->truncfactor = 0.0; 850*0f1074feSSatish Balay jac->agg_num_paths = 1; 85116d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr); 85216d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr); 85316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 85416d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 85516d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr); 85616d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr); 85716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr); 85816d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr); 85916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr); 86016d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr); 861*0f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);CHKERRQ(ierr); 862*0f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl); 863*0f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);CHKERRQ(ierr); 864*0f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);CHKERRQ(ierr); 865*0f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, jac->relaxtype[0]); /*defaults coarse to 9*/ 866*0f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver, jac->gridsweeps[0]); /*defaults coarse to 1 */ 867*0f1074feSSatish Balay 868*0f1074feSSatish Balay 869*0f1074feSSatish Balay 870*0f1074feSSatish Balay 87116d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 87216d9e3a6SLisandro Dalcin } 87316d9e3a6SLisandro Dalcin ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr); 87416d9e3a6SLisandro Dalcin jac->hypre_type = PETSC_NULL; 87516d9e3a6SLisandro Dalcin SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name); 87616d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 87716d9e3a6SLisandro Dalcin } 87816d9e3a6SLisandro Dalcin EXTERN_C_END 87916d9e3a6SLisandro Dalcin 88016d9e3a6SLisandro Dalcin /* 88116d9e3a6SLisandro Dalcin It only gets here if the HYPRE type has not been set before the call to 88216d9e3a6SLisandro Dalcin ...SetFromOptions() which actually is most of the time 88316d9e3a6SLisandro Dalcin */ 88416d9e3a6SLisandro Dalcin #undef __FUNCT__ 88516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE" 88616d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE(PC pc) 88716d9e3a6SLisandro Dalcin { 88816d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 88916d9e3a6SLisandro Dalcin PetscErrorCode ierr; 89016d9e3a6SLisandro Dalcin int indx; 89116d9e3a6SLisandro Dalcin const char *type[] = {"pilut","parasails","boomeramg","euclid"}; 89216d9e3a6SLisandro Dalcin PetscTruth flg; 89316d9e3a6SLisandro Dalcin 89416d9e3a6SLisandro Dalcin PetscFunctionBegin; 89516d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr); 89616d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"pilut",&indx,&flg);CHKERRQ(ierr); 89716d9e3a6SLisandro Dalcin if (PetscOptionsPublishCount) { /* force the default if it was not yet set and user did not set with option */ 89816d9e3a6SLisandro Dalcin if (!flg && !jac->hypre_type) { 89916d9e3a6SLisandro Dalcin flg = PETSC_TRUE; 90016d9e3a6SLisandro Dalcin indx = 0; 90116d9e3a6SLisandro Dalcin } 90216d9e3a6SLisandro Dalcin } 90316d9e3a6SLisandro Dalcin if (flg) { 90416d9e3a6SLisandro Dalcin ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr); 90516d9e3a6SLisandro Dalcin } 90616d9e3a6SLisandro Dalcin if (pc->ops->setfromoptions) { 90716d9e3a6SLisandro Dalcin ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr); 90816d9e3a6SLisandro Dalcin } 90916d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 91016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 91116d9e3a6SLisandro Dalcin } 91216d9e3a6SLisandro Dalcin 91316d9e3a6SLisandro Dalcin #undef __FUNCT__ 91416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType" 91516d9e3a6SLisandro Dalcin /*@C 91616d9e3a6SLisandro Dalcin PCHYPRESetType - Sets which hypre preconditioner you wish to use 91716d9e3a6SLisandro Dalcin 91816d9e3a6SLisandro Dalcin Input Parameters: 91916d9e3a6SLisandro Dalcin + pc - the preconditioner context 92016d9e3a6SLisandro Dalcin - name - either pilut, parasails, boomeramg, euclid 92116d9e3a6SLisandro Dalcin 92216d9e3a6SLisandro Dalcin Options Database Keys: 92316d9e3a6SLisandro Dalcin -pc_hypre_type - One of pilut, parasails, boomeramg, euclid 92416d9e3a6SLisandro Dalcin 92516d9e3a6SLisandro Dalcin Level: intermediate 92616d9e3a6SLisandro Dalcin 92716d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 92816d9e3a6SLisandro Dalcin PCHYPRE 92916d9e3a6SLisandro Dalcin 93016d9e3a6SLisandro Dalcin @*/ 93116d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC pc,const char name[]) 93216d9e3a6SLisandro Dalcin { 93316d9e3a6SLisandro Dalcin PetscErrorCode ierr,(*f)(PC,const char[]); 93416d9e3a6SLisandro Dalcin 93516d9e3a6SLisandro Dalcin PetscFunctionBegin; 93616d9e3a6SLisandro Dalcin PetscValidHeaderSpecific(pc,PC_COOKIE,1); 93716d9e3a6SLisandro Dalcin PetscValidCharPointer(name,2); 93816d9e3a6SLisandro Dalcin ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPRESetType_C",(void (**)(void))&f);CHKERRQ(ierr); 93916d9e3a6SLisandro Dalcin if (f) { 94016d9e3a6SLisandro Dalcin ierr = (*f)(pc,name);CHKERRQ(ierr); 94116d9e3a6SLisandro Dalcin } 94216d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 94316d9e3a6SLisandro Dalcin } 94416d9e3a6SLisandro Dalcin 94516d9e3a6SLisandro Dalcin #undef __FUNCT__ 94616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType" 94716d9e3a6SLisandro Dalcin /*@C 94816d9e3a6SLisandro Dalcin PCHYPREGetType - Gets which hypre preconditioner you are using 94916d9e3a6SLisandro Dalcin 95016d9e3a6SLisandro Dalcin Input Parameter: 95116d9e3a6SLisandro Dalcin . pc - the preconditioner context 95216d9e3a6SLisandro Dalcin 95316d9e3a6SLisandro Dalcin Output Parameter: 95416d9e3a6SLisandro Dalcin . name - either pilut, parasails, boomeramg, euclid 95516d9e3a6SLisandro Dalcin 95616d9e3a6SLisandro Dalcin Level: intermediate 95716d9e3a6SLisandro Dalcin 95816d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC, 95916d9e3a6SLisandro Dalcin PCHYPRE 96016d9e3a6SLisandro Dalcin 96116d9e3a6SLisandro Dalcin @*/ 96216d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC pc,const char *name[]) 96316d9e3a6SLisandro Dalcin { 96416d9e3a6SLisandro Dalcin PetscErrorCode ierr,(*f)(PC,const char*[]); 96516d9e3a6SLisandro Dalcin 96616d9e3a6SLisandro Dalcin PetscFunctionBegin; 96716d9e3a6SLisandro Dalcin PetscValidHeaderSpecific(pc,PC_COOKIE,1); 96816d9e3a6SLisandro Dalcin PetscValidPointer(name,2); 96916d9e3a6SLisandro Dalcin ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPREGetType_C",(void (**)(void))&f);CHKERRQ(ierr); 97016d9e3a6SLisandro Dalcin if (f) { 97116d9e3a6SLisandro Dalcin ierr = (*f)(pc,name);CHKERRQ(ierr); 97216d9e3a6SLisandro Dalcin } 97316d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 97416d9e3a6SLisandro Dalcin } 97516d9e3a6SLisandro Dalcin 97616d9e3a6SLisandro Dalcin /*MC 97716d9e3a6SLisandro Dalcin PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre 97816d9e3a6SLisandro Dalcin 97916d9e3a6SLisandro Dalcin Options Database Keys: 98016d9e3a6SLisandro Dalcin + -pc_hypre_type - One of pilut, parasails, boomeramg, euclid 98116d9e3a6SLisandro Dalcin - Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX 98216d9e3a6SLisandro Dalcin preconditioner 98316d9e3a6SLisandro Dalcin 98416d9e3a6SLisandro Dalcin Level: intermediate 98516d9e3a6SLisandro Dalcin 98616d9e3a6SLisandro Dalcin Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()), 98716d9e3a6SLisandro Dalcin the many hypre options can ONLY be set via the options database (e.g. the command line 98816d9e3a6SLisandro Dalcin or with PetscOptionsSetValue(), there are no functions to set them) 98916d9e3a6SLisandro Dalcin 99016d9e3a6SLisandro Dalcin The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations 991*0f1074feSSatish Balay (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if 992*0f1074feSSatish Balay -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner 993*0f1074feSSatish Balay (-pc_hypre_boomeramg_rtol should be set to 0.0 - the default - to strictly use a fixed number of 994*0f1074feSSatish Balay iterations per hypre call). -ksp_max_iter and -ksp_rtol STILL determine the total number of iterations 995*0f1074feSSatish Balay and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10 996*0f1074feSSatish Balay then AT MOST twenty V-cycles of boomeramg will be called. 99716d9e3a6SLisandro Dalcin 998*0f1074feSSatish Balay Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation 999*0f1074feSSatish Balay (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry. 1000*0f1074feSSatish Balay Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi. 100116d9e3a6SLisandro Dalcin If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly 100216d9e3a6SLisandro Dalcin and use -ksp_max_it to control the number of V-cycles. 100316d9e3a6SLisandro Dalcin (see the PETSc FAQ.html at the PETSc website under the Documentation tab). 100416d9e3a6SLisandro Dalcin 100516d9e3a6SLisandro Dalcin 2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option 100616d9e3a6SLisandro Dalcin -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L. 100716d9e3a6SLisandro Dalcin 100816d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 100916d9e3a6SLisandro Dalcin PCHYPRESetType() 101016d9e3a6SLisandro Dalcin 101116d9e3a6SLisandro Dalcin M*/ 101216d9e3a6SLisandro Dalcin 101316d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 101416d9e3a6SLisandro Dalcin #undef __FUNCT__ 101516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCCreate_HYPRE" 101616d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_HYPRE(PC pc) 101716d9e3a6SLisandro Dalcin { 101816d9e3a6SLisandro Dalcin PC_HYPRE *jac; 101916d9e3a6SLisandro Dalcin PetscErrorCode ierr; 102016d9e3a6SLisandro Dalcin 102116d9e3a6SLisandro Dalcin PetscFunctionBegin; 102216d9e3a6SLisandro Dalcin ierr = PetscNew(PC_HYPRE,&jac);CHKERRQ(ierr); 102316d9e3a6SLisandro Dalcin ierr = PetscLogObjectMemory(pc,sizeof(PC_HYPRE));CHKERRQ(ierr); 102416d9e3a6SLisandro Dalcin pc->data = jac; 102516d9e3a6SLisandro Dalcin pc->ops->destroy = PCDestroy_HYPRE; 102616d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE; 102716d9e3a6SLisandro Dalcin pc->ops->setup = PCSetUp_HYPRE; 102816d9e3a6SLisandro Dalcin pc->ops->apply = PCApply_HYPRE; 102916d9e3a6SLisandro Dalcin jac->comm_hypre = MPI_COMM_NULL; 103016d9e3a6SLisandro Dalcin jac->hypre_type = PETSC_NULL; 103116d9e3a6SLisandro Dalcin /* duplicate communicator for hypre */ 103216d9e3a6SLisandro Dalcin ierr = MPI_Comm_dup(pc->comm,&(jac->comm_hypre));CHKERRQ(ierr); 103316d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr); 103416d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr); 103516d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 103616d9e3a6SLisandro Dalcin } 103716d9e3a6SLisandro Dalcin EXTERN_C_END 1038