116d9e3a6SLisandro Dalcin #define PETSCKSP_DLL 216d9e3a6SLisandro Dalcin 316d9e3a6SLisandro Dalcin /* 416d9e3a6SLisandro Dalcin Provides an interface to the LLNL package hypre 516d9e3a6SLisandro Dalcin */ 60f1074feSSatish Balay 70f1074feSSatish Balay /* Must use hypre 2.0.0 or more recent. */ 80f1074feSSatish Balay 916d9e3a6SLisandro Dalcin #include "private/pcimpl.h" /*I "petscpc.h" I*/ 1016d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 1116d9e3a6SLisandro Dalcin #include "HYPRE.h" 120f1074feSSatish Balay #include "HYPRE_parcsr_ls.h" 130f1074feSSatish Balay #include "_hypre_parcsr_mv.h" 140f1074feSSatish 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; 650f1074feSSatish Balay int gridsweeps[3]; 6616d9e3a6SLisandro Dalcin int coarsentype; 6716d9e3a6SLisandro Dalcin int measuretype; 680f1074feSSatish Balay int relaxtype[3]; 6916d9e3a6SLisandro Dalcin double relaxweight; 7016d9e3a6SLisandro Dalcin double outerrelaxweight; 7116d9e3a6SLisandro Dalcin int relaxorder; 7216d9e3a6SLisandro Dalcin double truncfactor; 7316d9e3a6SLisandro Dalcin PetscTruth applyrichardson; 740f1074feSSatish Balay int pmax; 750f1074feSSatish Balay int interptype; 760f1074feSSatish Balay int agg_nl; 770f1074feSSatish Balay int agg_num_paths; 788f87f92bSBarry Smith int nodal_coarsen; 798f87f92bSBarry Smith PetscTruth nodal_relax; 808f87f92bSBarry Smith int nodal_relax_levels; 8116d9e3a6SLisandro Dalcin } PC_HYPRE; 8216d9e3a6SLisandro Dalcin 8316d9e3a6SLisandro Dalcin 8416d9e3a6SLisandro Dalcin #undef __FUNCT__ 8516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetUp_HYPRE" 8616d9e3a6SLisandro Dalcin static PetscErrorCode PCSetUp_HYPRE(PC pc) 8716d9e3a6SLisandro Dalcin { 8816d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 8916d9e3a6SLisandro Dalcin PetscErrorCode ierr; 9016d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 9116d9e3a6SLisandro Dalcin HYPRE_ParVector bv,xv; 9216d9e3a6SLisandro Dalcin PetscInt bs; 9316d9e3a6SLisandro Dalcin int hierr; 9416d9e3a6SLisandro Dalcin 9516d9e3a6SLisandro Dalcin PetscFunctionBegin; 9616d9e3a6SLisandro Dalcin if (!jac->hypre_type) { 9716d9e3a6SLisandro Dalcin ierr = PCHYPRESetType(pc,"pilut");CHKERRQ(ierr); 9816d9e3a6SLisandro Dalcin } 9916d9e3a6SLisandro Dalcin #if 1 10016d9e3a6SLisandro Dalcin if (!pc->setupcalled) { 10116d9e3a6SLisandro Dalcin /* create the matrix and vectors the first time through */ 10216d9e3a6SLisandro Dalcin Vec x,b; 10316d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr); 10416d9e3a6SLisandro Dalcin ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr); 10516d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr); 10616d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr); 10716d9e3a6SLisandro Dalcin ierr = VecDestroy(x);CHKERRQ(ierr); 10816d9e3a6SLisandro Dalcin ierr = VecDestroy(b);CHKERRQ(ierr); 10916d9e3a6SLisandro Dalcin } else if (pc->flag != SAME_NONZERO_PATTERN) { 11016d9e3a6SLisandro Dalcin /* rebuild the matrix from scratch */ 11116d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); 11216d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr); 11316d9e3a6SLisandro Dalcin } 11416d9e3a6SLisandro Dalcin #else 11516d9e3a6SLisandro Dalcin if (!jac->ij) { /* create the matrix the first time through */ 11616d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr); 11716d9e3a6SLisandro Dalcin } 11816d9e3a6SLisandro Dalcin if (!jac->b) { /* create the vectors the first time through */ 11916d9e3a6SLisandro Dalcin Vec x,b; 12016d9e3a6SLisandro Dalcin ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr); 12116d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr); 12216d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr); 12316d9e3a6SLisandro Dalcin ierr = VecDestroy(x);CHKERRQ(ierr); 12416d9e3a6SLisandro Dalcin ierr = VecDestroy(b);CHKERRQ(ierr); 12516d9e3a6SLisandro Dalcin } 12616d9e3a6SLisandro Dalcin #endif 12716d9e3a6SLisandro Dalcin /* special case for BoomerAMG */ 12816d9e3a6SLisandro Dalcin if (jac->setup == HYPRE_BoomerAMGSetup) { 12916d9e3a6SLisandro Dalcin ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr); 13016d9e3a6SLisandro Dalcin if (bs > 1) { 13116d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetNumFunctions(jac->hsolver,bs);CHKERRQ(ierr); 13216d9e3a6SLisandro Dalcin } 13316d9e3a6SLisandro Dalcin }; 13416d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr); 13516d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 13616d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&bv);CHKERRQ(ierr); 13716d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&xv);CHKERRQ(ierr); 13816d9e3a6SLisandro Dalcin hierr = (*jac->setup)(jac->hsolver,hmat,bv,xv); 13916d9e3a6SLisandro Dalcin if (hierr) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE setup, error code %d",hierr); 14016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 14116d9e3a6SLisandro Dalcin } 14216d9e3a6SLisandro Dalcin 14316d9e3a6SLisandro Dalcin /* 14416d9e3a6SLisandro Dalcin Replaces the address where the HYPRE vector points to its data with the address of 14516d9e3a6SLisandro Dalcin PETSc's data. Saves the old address so it can be reset when we are finished with it. 14616d9e3a6SLisandro Dalcin Allows use to get the data into a HYPRE vector without the cost of memcopies 14716d9e3a6SLisandro Dalcin */ 14816d9e3a6SLisandro Dalcin #define HYPREReplacePointer(b,newvalue,savedvalue) {\ 14916d9e3a6SLisandro Dalcin hypre_ParVector *par_vector = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));\ 15016d9e3a6SLisandro Dalcin hypre_Vector *local_vector = hypre_ParVectorLocalVector(par_vector);\ 15116d9e3a6SLisandro Dalcin savedvalue = local_vector->data;\ 15216d9e3a6SLisandro Dalcin local_vector->data = newvalue;} 15316d9e3a6SLisandro Dalcin 15416d9e3a6SLisandro Dalcin #undef __FUNCT__ 15516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApply_HYPRE" 15616d9e3a6SLisandro Dalcin static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x) 15716d9e3a6SLisandro Dalcin { 15816d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 15916d9e3a6SLisandro Dalcin PetscErrorCode ierr; 16016d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 16116d9e3a6SLisandro Dalcin PetscScalar *bv,*xv; 16216d9e3a6SLisandro Dalcin HYPRE_ParVector jbv,jxv; 16316d9e3a6SLisandro Dalcin PetscScalar *sbv,*sxv; 16416d9e3a6SLisandro Dalcin int hierr; 16516d9e3a6SLisandro Dalcin 16616d9e3a6SLisandro Dalcin PetscFunctionBegin; 16716d9e3a6SLisandro Dalcin if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);} 16816d9e3a6SLisandro Dalcin ierr = VecGetArray(b,&bv);CHKERRQ(ierr); 16916d9e3a6SLisandro Dalcin ierr = VecGetArray(x,&xv);CHKERRQ(ierr); 17016d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,bv,sbv); 17116d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,xv,sxv); 17216d9e3a6SLisandro Dalcin 17316d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 17416d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr); 17516d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr); 17616d9e3a6SLisandro Dalcin hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv); 1770f1074feSSatish Balay 1780f1074feSSatish Balay /*if (hierr && (hierr != HYPRE_ERROR_CONV || jac->solve != HYPRE_BoomerAMGSolve))SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 1790f1074feSSatish Balay */ 1800f1074feSSatish Balay /* error code of HYPRE_ERROR_CONV means convergence not achieved - if 1810f1074feSSatish Balay the tolerance is set to 0.0 (the default), a convergence error will 1820f1074feSSatish Balay not occur (so we may not want to overide the conv. error here?*/ 1830f1074feSSatish Balay if (hierr && hierr != HYPRE_ERROR_CONV) 1840f1074feSSatish Balay { 1850f1074feSSatish Balay SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 1860f1074feSSatish Balay } 18716d9e3a6SLisandro Dalcin if (hierr) hypre__global_error = 0; 18816d9e3a6SLisandro Dalcin 1890f1074feSSatish Balay 19016d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,sbv,bv); 19116d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,sxv,xv); 19216d9e3a6SLisandro Dalcin ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr); 19316d9e3a6SLisandro Dalcin ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr); 19416d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 19516d9e3a6SLisandro Dalcin } 19616d9e3a6SLisandro Dalcin 19716d9e3a6SLisandro Dalcin #undef __FUNCT__ 19816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCDestroy_HYPRE" 19916d9e3a6SLisandro Dalcin static PetscErrorCode PCDestroy_HYPRE(PC pc) 20016d9e3a6SLisandro Dalcin { 20116d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 20216d9e3a6SLisandro Dalcin PetscErrorCode ierr; 20316d9e3a6SLisandro Dalcin 20416d9e3a6SLisandro Dalcin PetscFunctionBegin; 20516d9e3a6SLisandro Dalcin if (jac->ij) { ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); } 20616d9e3a6SLisandro Dalcin if (jac->b) { ierr = HYPRE_IJVectorDestroy(jac->b);CHKERRQ(ierr); } 20716d9e3a6SLisandro Dalcin if (jac->x) { ierr = HYPRE_IJVectorDestroy(jac->x);CHKERRQ(ierr); } 20816d9e3a6SLisandro Dalcin if (jac->destroy) { ierr = (*jac->destroy)(jac->hsolver);CHKERRQ(ierr); } 20916d9e3a6SLisandro Dalcin ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr); 21016d9e3a6SLisandro Dalcin if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);} 21116d9e3a6SLisandro Dalcin ierr = PetscFree(jac);CHKERRQ(ierr); 21216d9e3a6SLisandro Dalcin 21316d9e3a6SLisandro Dalcin ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr); 21416d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);CHKERRQ(ierr); 21516d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);CHKERRQ(ierr); 21616d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 21716d9e3a6SLisandro Dalcin } 21816d9e3a6SLisandro Dalcin 21916d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 22016d9e3a6SLisandro Dalcin #undef __FUNCT__ 22116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut" 22216d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc) 22316d9e3a6SLisandro Dalcin { 22416d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 22516d9e3a6SLisandro Dalcin PetscErrorCode ierr; 22616d9e3a6SLisandro Dalcin PetscTruth flag; 22716d9e3a6SLisandro Dalcin 22816d9e3a6SLisandro Dalcin PetscFunctionBegin; 22916d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr); 23016d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr); 23116d9e3a6SLisandro Dalcin if (flag) { 23216d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 23316d9e3a6SLisandro Dalcin } 23416d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr); 23516d9e3a6SLisandro Dalcin if (flag) { 23616d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutSetDropTolerance(jac->hsolver,jac->tol);CHKERRQ(ierr); 23716d9e3a6SLisandro Dalcin } 23816d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr); 23916d9e3a6SLisandro Dalcin if (flag) { 24016d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutSetFactorRowSize(jac->hsolver,jac->factorrowsize);CHKERRQ(ierr); 24116d9e3a6SLisandro Dalcin } 24216d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 24316d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 24416d9e3a6SLisandro Dalcin } 24516d9e3a6SLisandro Dalcin 24616d9e3a6SLisandro Dalcin #undef __FUNCT__ 24716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Pilut" 24816d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer) 24916d9e3a6SLisandro Dalcin { 25016d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 25116d9e3a6SLisandro Dalcin PetscErrorCode ierr; 25216d9e3a6SLisandro Dalcin PetscTruth iascii; 25316d9e3a6SLisandro Dalcin 25416d9e3a6SLisandro Dalcin PetscFunctionBegin; 25516d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 25616d9e3a6SLisandro Dalcin if (iascii) { 25716d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut preconditioning\n");CHKERRQ(ierr); 25816d9e3a6SLisandro Dalcin if (jac->maxiter != PETSC_DEFAULT) { 25916d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr); 26016d9e3a6SLisandro Dalcin } else { 26116d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr); 26216d9e3a6SLisandro Dalcin } 26316d9e3a6SLisandro Dalcin if (jac->tol != PETSC_DEFAULT) { 26416d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr); 26516d9e3a6SLisandro Dalcin } else { 26616d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr); 26716d9e3a6SLisandro Dalcin } 26816d9e3a6SLisandro Dalcin if (jac->factorrowsize != PETSC_DEFAULT) { 26916d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr); 27016d9e3a6SLisandro Dalcin } else { 27116d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default factor row size \n");CHKERRQ(ierr); 27216d9e3a6SLisandro Dalcin } 27316d9e3a6SLisandro Dalcin } 27416d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 27516d9e3a6SLisandro Dalcin } 27616d9e3a6SLisandro Dalcin 27716d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 27816d9e3a6SLisandro Dalcin #undef __FUNCT__ 27916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid" 28016d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc) 28116d9e3a6SLisandro Dalcin { 28216d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 28316d9e3a6SLisandro Dalcin PetscErrorCode ierr; 28416d9e3a6SLisandro Dalcin PetscTruth flag; 285390e7148SBarry Smith char *args[8],levels[16]; 286390e7148SBarry Smith PetscInt cnt = 0; 28716d9e3a6SLisandro Dalcin 28816d9e3a6SLisandro Dalcin PetscFunctionBegin; 28916d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr); 29016d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr); 29116d9e3a6SLisandro Dalcin if (flag) { 29216d9e3a6SLisandro Dalcin if (jac->levels < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels); 29316d9e3a6SLisandro Dalcin sprintf(levels,"%d",jac->levels); 294390e7148SBarry Smith args[cnt++] = (char*)"-level"; args[cnt++] = levels; 29516d9e3a6SLisandro Dalcin } 29616d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);CHKERRQ(ierr); 29716d9e3a6SLisandro Dalcin if (jac->bjilu) { 298390e7148SBarry Smith args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1"; 29916d9e3a6SLisandro Dalcin } 30016d9e3a6SLisandro Dalcin 30116d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);CHKERRQ(ierr); 30216d9e3a6SLisandro Dalcin if (jac->printstatistics) { 303390e7148SBarry Smith args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1"; 304390e7148SBarry Smith args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1"; 30516d9e3a6SLisandro Dalcin } 30616d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 307390e7148SBarry Smith if (cnt) { 308390e7148SBarry Smith ierr = HYPRE_EuclidSetParams(jac->hsolver,cnt,args);CHKERRQ(ierr); 309390e7148SBarry Smith } 31016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 31116d9e3a6SLisandro Dalcin } 31216d9e3a6SLisandro Dalcin 31316d9e3a6SLisandro Dalcin #undef __FUNCT__ 31416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Euclid" 31516d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer) 31616d9e3a6SLisandro Dalcin { 31716d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 31816d9e3a6SLisandro Dalcin PetscErrorCode ierr; 31916d9e3a6SLisandro Dalcin PetscTruth iascii; 32016d9e3a6SLisandro Dalcin 32116d9e3a6SLisandro Dalcin PetscFunctionBegin; 32216d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 32316d9e3a6SLisandro Dalcin if (iascii) { 32416d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid preconditioning\n");CHKERRQ(ierr); 32516d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr); 32616d9e3a6SLisandro Dalcin if (jac->bjilu) { 32716d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr); 32816d9e3a6SLisandro Dalcin } 32916d9e3a6SLisandro Dalcin } 33016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 33116d9e3a6SLisandro Dalcin } 33216d9e3a6SLisandro Dalcin 33316d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 33416d9e3a6SLisandro Dalcin 33516d9e3a6SLisandro Dalcin #undef __FUNCT__ 33616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG" 33716d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x) 33816d9e3a6SLisandro Dalcin { 33916d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 34016d9e3a6SLisandro Dalcin PetscErrorCode ierr; 34116d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 34216d9e3a6SLisandro Dalcin PetscScalar *bv,*xv; 34316d9e3a6SLisandro Dalcin HYPRE_ParVector jbv,jxv; 34416d9e3a6SLisandro Dalcin PetscScalar *sbv,*sxv; 34516d9e3a6SLisandro Dalcin int hierr; 34616d9e3a6SLisandro Dalcin 34716d9e3a6SLisandro Dalcin PetscFunctionBegin; 34816d9e3a6SLisandro Dalcin ierr = VecSet(x,0.0);CHKERRQ(ierr); 34916d9e3a6SLisandro Dalcin ierr = VecGetArray(b,&bv);CHKERRQ(ierr); 35016d9e3a6SLisandro Dalcin ierr = VecGetArray(x,&xv);CHKERRQ(ierr); 35116d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,bv,sbv); 35216d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,xv,sxv); 35316d9e3a6SLisandro Dalcin 35416d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 35516d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr); 35616d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr); 35716d9e3a6SLisandro Dalcin 35816d9e3a6SLisandro Dalcin hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv); 35916d9e3a6SLisandro Dalcin /* error code of 1 in BoomerAMG merely means convergence not achieved */ 36016d9e3a6SLisandro Dalcin if (hierr && (hierr != 1)) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 36116d9e3a6SLisandro Dalcin if (hierr) hypre__global_error = 0; 36216d9e3a6SLisandro Dalcin 36316d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,sbv,bv); 36416d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,sxv,xv); 36516d9e3a6SLisandro Dalcin ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr); 36616d9e3a6SLisandro Dalcin ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr); 36716d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 36816d9e3a6SLisandro Dalcin } 36916d9e3a6SLisandro Dalcin 37016d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCycleType[] = {"","V","W"}; 3710f1074feSSatish Balay static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"}; 37216d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGMeasureType[] = {"local","global"}; 3730f1074feSSatish Balay static const char *HYPREBoomerAMGRelaxType[] = {"Jacobi","sequential-Gauss-Seidel","","SOR/Jacobi","backward-SOR/Jacobi","","symmetric-SOR/Jacobi", 3740f1074feSSatish Balay "","","Gaussian-elimination"}; 3750f1074feSSatish Balay static const char *HYPREBoomerAMGInterpType[] = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i", 3760f1074feSSatish Balay "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"}; 37716d9e3a6SLisandro Dalcin #undef __FUNCT__ 37816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG" 37916d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc) 38016d9e3a6SLisandro Dalcin { 38116d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 38216d9e3a6SLisandro Dalcin PetscErrorCode ierr; 38316d9e3a6SLisandro Dalcin int n,indx; 38416d9e3a6SLisandro Dalcin PetscTruth flg, tmp_truth; 38516d9e3a6SLisandro Dalcin double tmpdbl, twodbl[2]; 38616d9e3a6SLisandro Dalcin 38716d9e3a6SLisandro Dalcin PetscFunctionBegin; 38816d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr); 38916d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr); 39016d9e3a6SLisandro Dalcin if (flg) { 39116d9e3a6SLisandro Dalcin jac->cycletype = indx; 39216d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr); 39316d9e3a6SLisandro Dalcin } 39416d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr); 39516d9e3a6SLisandro Dalcin if (flg) { 39616d9e3a6SLisandro Dalcin if (jac->maxlevels < 2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels); 39716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr); 39816d9e3a6SLisandro Dalcin } 39916d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr); 40016d9e3a6SLisandro Dalcin if (flg) { 40116d9e3a6SLisandro Dalcin if (jac->maxiter < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter); 40216d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 40316d9e3a6SLisandro Dalcin } 4040f1074feSSatish 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); 40516d9e3a6SLisandro Dalcin if (flg) { 4060f1074feSSatish Balay if (jac->tol < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be greater than or equal to zero",jac->tol); 40716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 40816d9e3a6SLisandro Dalcin } 40916d9e3a6SLisandro Dalcin 4100f1074feSSatish Balay ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr); 41116d9e3a6SLisandro Dalcin if (flg) { 41216d9e3a6SLisandro Dalcin if (jac->truncfactor < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor); 41316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr); 41416d9e3a6SLisandro Dalcin } 41516d9e3a6SLisandro Dalcin 4160f1074feSSatish 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); 4170f1074feSSatish Balay if (flg) { 4180f1074feSSatish Balay if (jac->pmax < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"P_max %G must be greater than or equal to zero",jac->pmax); 4190f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);CHKERRQ(ierr); 4200f1074feSSatish Balay } 4210f1074feSSatish Balay 4220f1074feSSatish Balay ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);CHKERRQ(ierr); 4230f1074feSSatish Balay if (flg) { 4240f1074feSSatish 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); 4250f1074feSSatish Balay 4260f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl);CHKERRQ(ierr); 4270f1074feSSatish Balay } 4280f1074feSSatish Balay 4290f1074feSSatish Balay 4300f1074feSSatish 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); 4310f1074feSSatish Balay if (flg) { 4320f1074feSSatish 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); 4330f1074feSSatish Balay 4340f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);CHKERRQ(ierr); 4350f1074feSSatish Balay } 4360f1074feSSatish Balay 4370f1074feSSatish Balay 43816d9e3a6SLisandro Dalcin ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr); 43916d9e3a6SLisandro Dalcin if (flg) { 44016d9e3a6SLisandro Dalcin if (jac->strongthreshold < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold); 44116d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr); 44216d9e3a6SLisandro Dalcin } 44316d9e3a6SLisandro Dalcin ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr); 44416d9e3a6SLisandro Dalcin if (flg) { 44516d9e3a6SLisandro Dalcin if (jac->maxrowsum < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum); 44616d9e3a6SLisandro Dalcin if (jac->maxrowsum > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum); 44716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr); 44816d9e3a6SLisandro Dalcin } 44916d9e3a6SLisandro Dalcin 45016d9e3a6SLisandro Dalcin /* Grid sweeps */ 4510f1074feSSatish 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); 45216d9e3a6SLisandro Dalcin if (flg) { 45316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver,indx);CHKERRQ(ierr); 45416d9e3a6SLisandro Dalcin /* modify the jac structure so we can view the updated options with PC_View */ 45516d9e3a6SLisandro Dalcin jac->gridsweeps[0] = indx; 4560f1074feSSatish Balay jac->gridsweeps[1] = indx; 4570f1074feSSatish Balay /*defaults coarse to 1 */ 4580f1074feSSatish Balay jac->gridsweeps[2] = 1; 45916d9e3a6SLisandro Dalcin } 4600f1074feSSatish Balay 4610f1074feSSatish Balay ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx ,&flg);CHKERRQ(ierr); 46216d9e3a6SLisandro Dalcin if (flg) { 46316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 1);CHKERRQ(ierr); 4640f1074feSSatish Balay jac->gridsweeps[0] = indx; 46516d9e3a6SLisandro Dalcin } 46616d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr); 46716d9e3a6SLisandro Dalcin if (flg) { 46816d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 2);CHKERRQ(ierr); 4690f1074feSSatish Balay jac->gridsweeps[1] = indx; 47016d9e3a6SLisandro Dalcin } 4710f1074feSSatish Balay ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr); 47216d9e3a6SLisandro Dalcin if (flg) { 47316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 3);CHKERRQ(ierr); 4740f1074feSSatish Balay jac->gridsweeps[2] = indx; 47516d9e3a6SLisandro Dalcin } 47616d9e3a6SLisandro Dalcin 47716d9e3a6SLisandro Dalcin /* Relax type */ 4780f1074feSSatish 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); 47916d9e3a6SLisandro Dalcin if (flg) { 4800f1074feSSatish Balay jac->relaxtype[0] = jac->relaxtype[1] = indx; 48116d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, indx);CHKERRQ(ierr); 4820f1074feSSatish Balay /* by default, coarse type set to 9 */ 4830f1074feSSatish Balay jac->relaxtype[2] = 9; 4840f1074feSSatish Balay 48516d9e3a6SLisandro Dalcin } 4860f1074feSSatish Balay ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr); 48716d9e3a6SLisandro Dalcin if (flg) { 48816d9e3a6SLisandro Dalcin jac->relaxtype[0] = indx; 48916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 1);CHKERRQ(ierr); 49016d9e3a6SLisandro Dalcin } 4910f1074feSSatish Balay ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr); 49216d9e3a6SLisandro Dalcin if (flg) { 4930f1074feSSatish Balay jac->relaxtype[1] = indx; 49416d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 2);CHKERRQ(ierr); 49516d9e3a6SLisandro Dalcin } 49616d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr); 49716d9e3a6SLisandro Dalcin if (flg) { 4980f1074feSSatish Balay jac->relaxtype[2] = indx; 49916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 3);CHKERRQ(ierr); 50016d9e3a6SLisandro Dalcin } 50116d9e3a6SLisandro Dalcin 50216d9e3a6SLisandro Dalcin /* Relaxation Weight */ 50316d9e3a6SLisandro 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); 50416d9e3a6SLisandro Dalcin if (flg) { 50516d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxWt(jac->hsolver,tmpdbl);CHKERRQ(ierr); 50616d9e3a6SLisandro Dalcin jac->relaxweight = tmpdbl; 50716d9e3a6SLisandro Dalcin } 50816d9e3a6SLisandro Dalcin 50916d9e3a6SLisandro Dalcin n=2; 51016d9e3a6SLisandro Dalcin twodbl[0] = twodbl[1] = 1.0; 51116d9e3a6SLisandro 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); 51216d9e3a6SLisandro Dalcin if (flg) { 51316d9e3a6SLisandro Dalcin if (n == 2) { 51416d9e3a6SLisandro Dalcin indx = (int)PetscAbsReal(twodbl[1]); 51516d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetLevelRelaxWt(jac->hsolver,twodbl[0],indx);CHKERRQ(ierr); 51616d9e3a6SLisandro Dalcin } else { 51716d9e3a6SLisandro 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); 51816d9e3a6SLisandro Dalcin } 51916d9e3a6SLisandro Dalcin } 52016d9e3a6SLisandro Dalcin 52116d9e3a6SLisandro Dalcin /* Outer relaxation Weight */ 52216d9e3a6SLisandro 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); 52316d9e3a6SLisandro Dalcin if (flg) { 52416d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetOuterWt( jac->hsolver, tmpdbl);CHKERRQ(ierr); 52516d9e3a6SLisandro Dalcin jac->outerrelaxweight = tmpdbl; 52616d9e3a6SLisandro Dalcin } 52716d9e3a6SLisandro Dalcin 52816d9e3a6SLisandro Dalcin n=2; 52916d9e3a6SLisandro Dalcin twodbl[0] = twodbl[1] = 1.0; 53016d9e3a6SLisandro 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); 53116d9e3a6SLisandro Dalcin if (flg) { 53216d9e3a6SLisandro Dalcin if (n == 2) { 53316d9e3a6SLisandro Dalcin indx = (int)PetscAbsReal(twodbl[1]); 53416d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetLevelOuterWt( jac->hsolver, twodbl[0], indx);CHKERRQ(ierr); 53516d9e3a6SLisandro Dalcin } else { 53616d9e3a6SLisandro 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); 53716d9e3a6SLisandro Dalcin } 53816d9e3a6SLisandro Dalcin } 53916d9e3a6SLisandro Dalcin 54016d9e3a6SLisandro Dalcin /* the Relax Order */ 54116d9e3a6SLisandro Dalcin /* ierr = PetscOptionsName("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", &flg);CHKERRQ(ierr); */ 54216d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr); 54316d9e3a6SLisandro Dalcin 54416d9e3a6SLisandro Dalcin if (flg) { 54516d9e3a6SLisandro Dalcin jac->relaxorder = 0; 54616d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr); 54716d9e3a6SLisandro Dalcin } 54816d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr); 54916d9e3a6SLisandro Dalcin if (flg) { 55016d9e3a6SLisandro Dalcin jac->measuretype = indx; 55116d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr); 55216d9e3a6SLisandro Dalcin } 5530f1074feSSatish Balay /* update list length 3/07 */ 5540f1074feSSatish Balay ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,11,HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr); 55516d9e3a6SLisandro Dalcin if (flg) { 55616d9e3a6SLisandro Dalcin jac->coarsentype = indx; 55716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr); 55816d9e3a6SLisandro Dalcin } 5590f1074feSSatish Balay 5600f1074feSSatish Balay /* new 3/07 */ 5610f1074feSSatish Balay ierr = PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,14,HYPREBoomerAMGInterpType[0],&indx,&flg);CHKERRQ(ierr); 5620f1074feSSatish Balay if (flg) { 5630f1074feSSatish Balay jac->interptype = indx; 5640f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);CHKERRQ(ierr); 5650f1074feSSatish Balay } 5660f1074feSSatish Balay 5670f1074feSSatish Balay 5680f1074feSSatish Balay 56916d9e3a6SLisandro Dalcin ierr = PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);CHKERRQ(ierr); 57016d9e3a6SLisandro Dalcin if (flg) { 57116d9e3a6SLisandro Dalcin int level=3; 57216d9e3a6SLisandro Dalcin jac->printstatistics = PETSC_TRUE; 57316d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);CHKERRQ(ierr); 57416d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetPrintLevel(jac->hsolver,level);CHKERRQ(ierr); 57516d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetDebugFlag(jac->hsolver,level);CHKERRQ(ierr); 57616d9e3a6SLisandro Dalcin } 5778f87f92bSBarry Smith 5788f87f92bSBarry Smith ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr); 5798f87f92bSBarry Smith if (flg && tmp_truth) { 5808f87f92bSBarry Smith jac->nodal_coarsen = 1; 5818f87f92bSBarry Smith ierr = HYPRE_BoomerAMGSetNodal(jac->hsolver,1);CHKERRQ(ierr); 5828f87f92bSBarry Smith } 5838f87f92bSBarry Smith 5848f87f92bSBarry Smith ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr); 5858f87f92bSBarry Smith if (flg && tmp_truth) { 5868f87f92bSBarry Smith PetscInt tmp_int; 5878f87f92bSBarry Smith ierr = PetscOptionsInt( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);CHKERRQ(ierr); 5888f87f92bSBarry Smith if (flg) jac->nodal_relax_levels = tmp_int; 5898f87f92bSBarry Smith ierr = HYPRE_BoomerAMGSetSmoothType(jac->hsolver,6);CHKERRQ(ierr); 5908f87f92bSBarry Smith ierr = HYPRE_BoomerAMGSetDomainType(jac->hsolver,1);CHKERRQ(ierr); 5918f87f92bSBarry Smith ierr = HYPRE_BoomerAMGSetOverlap(jac->hsolver,0);CHKERRQ(ierr); 5928f87f92bSBarry Smith ierr = HYPRE_BoomerAMGSetSmoothNumLevels(jac->hsolver,jac->nodal_relax_levels);CHKERRQ(ierr); 5938f87f92bSBarry Smith } 5948f87f92bSBarry Smith 59516d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 59616d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 59716d9e3a6SLisandro Dalcin } 59816d9e3a6SLisandro Dalcin 59916d9e3a6SLisandro Dalcin #undef __FUNCT__ 60016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG" 60116d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its) 60216d9e3a6SLisandro Dalcin { 60316d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 60416d9e3a6SLisandro Dalcin PetscErrorCode ierr; 60516d9e3a6SLisandro Dalcin 60616d9e3a6SLisandro Dalcin PetscFunctionBegin; 60716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,its*jac->maxiter);CHKERRQ(ierr); 60816d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,PetscMin(rtol,jac->tol));CHKERRQ(ierr); 60916d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_TRUE; 61016d9e3a6SLisandro Dalcin ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr); 61116d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_FALSE; 61216d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 61316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 61416d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 61516d9e3a6SLisandro Dalcin } 61616d9e3a6SLisandro Dalcin 61716d9e3a6SLisandro Dalcin 61816d9e3a6SLisandro Dalcin #undef __FUNCT__ 61916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_BoomerAMG" 62016d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer) 62116d9e3a6SLisandro Dalcin { 62216d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 62316d9e3a6SLisandro Dalcin PetscErrorCode ierr; 62416d9e3a6SLisandro Dalcin PetscTruth iascii; 62516d9e3a6SLisandro Dalcin 62616d9e3a6SLisandro Dalcin PetscFunctionBegin; 62716d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 62816d9e3a6SLisandro Dalcin if (iascii) { 62916d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr); 63016d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr); 63116d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr); 63216d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr); 63316d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr); 63416d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr); 6350f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation truncation factor %G\n",jac->truncfactor);CHKERRQ(ierr); 6360f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);CHKERRQ(ierr); 6370f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);CHKERRQ(ierr); 6380f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);CHKERRQ(ierr); 6390f1074feSSatish Balay 64016d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr); 64116d9e3a6SLisandro Dalcin 6420f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps down %d\n",jac->gridsweeps[0]);CHKERRQ(ierr); 6430f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps up %d\n",jac->gridsweeps[1]);CHKERRQ(ierr); 6440f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps on coarse %d\n",jac->gridsweeps[2]);CHKERRQ(ierr); 64516d9e3a6SLisandro Dalcin 6460f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax down %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr); 6470f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax up %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr); 6480f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax on coarse %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr); 64916d9e3a6SLisandro Dalcin 65016d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax weight (all) %G\n",jac->relaxweight);CHKERRQ(ierr); 65116d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr); 65216d9e3a6SLisandro Dalcin 65316d9e3a6SLisandro Dalcin if (jac->relaxorder) { 65416d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr); 65516d9e3a6SLisandro Dalcin } else { 65616d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr); 65716d9e3a6SLisandro Dalcin } 65816d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Measure type %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr); 65916d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Coarsen type %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr); 6600f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation type %s\n",HYPREBoomerAMGInterpType[jac->interptype]);CHKERRQ(ierr); 6618f87f92bSBarry Smith if (jac->nodal_coarsen) { 6628f87f92bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal())\n");CHKERRQ(ierr); 6638f87f92bSBarry Smith } 6648f87f92bSBarry Smith if (jac->nodal_relax) { 6658f87f92bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);CHKERRQ(ierr); 6668f87f92bSBarry Smith } 66716d9e3a6SLisandro Dalcin } 66816d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 66916d9e3a6SLisandro Dalcin } 67016d9e3a6SLisandro Dalcin 67116d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 67216d9e3a6SLisandro Dalcin #undef __FUNCT__ 67316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails" 67416d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc) 67516d9e3a6SLisandro Dalcin { 67616d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 67716d9e3a6SLisandro Dalcin PetscErrorCode ierr; 67816d9e3a6SLisandro Dalcin int indx; 67916d9e3a6SLisandro Dalcin PetscTruth flag; 68016d9e3a6SLisandro Dalcin const char *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"}; 68116d9e3a6SLisandro Dalcin 68216d9e3a6SLisandro Dalcin PetscFunctionBegin; 68316d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr); 68416d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr); 68516d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr); 68616d9e3a6SLisandro Dalcin if (flag) { 68716d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr); 68816d9e3a6SLisandro Dalcin } 68916d9e3a6SLisandro Dalcin 69016d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr); 69116d9e3a6SLisandro Dalcin if (flag) { 69216d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr); 69316d9e3a6SLisandro Dalcin } 69416d9e3a6SLisandro Dalcin 69516d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr); 69616d9e3a6SLisandro Dalcin if (flag) { 69716d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr); 69816d9e3a6SLisandro Dalcin } 69916d9e3a6SLisandro Dalcin 70016d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,&flag);CHKERRQ(ierr); 70116d9e3a6SLisandro Dalcin if (flag) { 70216d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr); 70316d9e3a6SLisandro Dalcin } 70416d9e3a6SLisandro Dalcin 70516d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,&flag);CHKERRQ(ierr); 70616d9e3a6SLisandro Dalcin if (flag) { 70716d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr); 70816d9e3a6SLisandro Dalcin } 70916d9e3a6SLisandro Dalcin 71016d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);CHKERRQ(ierr); 71116d9e3a6SLisandro Dalcin if (flag) { 71216d9e3a6SLisandro Dalcin jac->symt = indx; 71316d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr); 71416d9e3a6SLisandro Dalcin } 71516d9e3a6SLisandro Dalcin 71616d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 71716d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 71816d9e3a6SLisandro Dalcin } 71916d9e3a6SLisandro Dalcin 72016d9e3a6SLisandro Dalcin #undef __FUNCT__ 72116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_ParaSails" 72216d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer) 72316d9e3a6SLisandro Dalcin { 72416d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 72516d9e3a6SLisandro Dalcin PetscErrorCode ierr; 72616d9e3a6SLisandro Dalcin PetscTruth iascii; 72716d9e3a6SLisandro Dalcin const char *symt = 0;; 72816d9e3a6SLisandro Dalcin 72916d9e3a6SLisandro Dalcin PetscFunctionBegin; 73016d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 73116d9e3a6SLisandro Dalcin if (iascii) { 73216d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails preconditioning\n");CHKERRQ(ierr); 73316d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr); 73416d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr); 73516d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr); 73616d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr); 73716d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: reuse nonzero structure %s\n",PetscTruths[jac->ruse]);CHKERRQ(ierr); 73816d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: print info to screen %s\n",PetscTruths[jac->logging]);CHKERRQ(ierr); 73916d9e3a6SLisandro Dalcin if (!jac->symt) { 74016d9e3a6SLisandro Dalcin symt = "nonsymmetric matrix and preconditioner"; 74116d9e3a6SLisandro Dalcin } else if (jac->symt == 1) { 74216d9e3a6SLisandro Dalcin symt = "SPD matrix and preconditioner"; 74316d9e3a6SLisandro Dalcin } else if (jac->symt == 2) { 74416d9e3a6SLisandro Dalcin symt = "nonsymmetric matrix but SPD preconditioner"; 74516d9e3a6SLisandro Dalcin } else { 74616d9e3a6SLisandro Dalcin SETERRQ1(PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt); 74716d9e3a6SLisandro Dalcin } 74816d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr); 74916d9e3a6SLisandro Dalcin } 75016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 75116d9e3a6SLisandro Dalcin } 75216d9e3a6SLisandro Dalcin /* ---------------------------------------------------------------------------------*/ 75316d9e3a6SLisandro Dalcin 75416d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 75516d9e3a6SLisandro Dalcin #undef __FUNCT__ 75616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType_HYPRE" 75716d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType_HYPRE(PC pc,const char *name[]) 75816d9e3a6SLisandro Dalcin { 75916d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 76016d9e3a6SLisandro Dalcin 76116d9e3a6SLisandro Dalcin PetscFunctionBegin; 76216d9e3a6SLisandro Dalcin *name = jac->hypre_type; 76316d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 76416d9e3a6SLisandro Dalcin } 76516d9e3a6SLisandro Dalcin EXTERN_C_END 76616d9e3a6SLisandro Dalcin 76716d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 76816d9e3a6SLisandro Dalcin #undef __FUNCT__ 76916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType_HYPRE" 77016d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType_HYPRE(PC pc,const char name[]) 77116d9e3a6SLisandro Dalcin { 77216d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 77316d9e3a6SLisandro Dalcin PetscErrorCode ierr; 77416d9e3a6SLisandro Dalcin PetscTruth flag; 77516d9e3a6SLisandro Dalcin 77616d9e3a6SLisandro Dalcin PetscFunctionBegin; 77716d9e3a6SLisandro Dalcin if (jac->hypre_type) { 77816d9e3a6SLisandro Dalcin ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr); 77916d9e3a6SLisandro Dalcin if (!flag) SETERRQ(PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set"); 78016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 78116d9e3a6SLisandro Dalcin } else { 78216d9e3a6SLisandro Dalcin ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr); 78316d9e3a6SLisandro Dalcin } 78416d9e3a6SLisandro Dalcin 78516d9e3a6SLisandro Dalcin jac->maxiter = PETSC_DEFAULT; 78616d9e3a6SLisandro Dalcin jac->tol = PETSC_DEFAULT; 78716d9e3a6SLisandro Dalcin jac->printstatistics = PetscLogPrintInfo; 78816d9e3a6SLisandro Dalcin 78916d9e3a6SLisandro Dalcin ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr); 79016d9e3a6SLisandro Dalcin if (flag) { 79116d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutCreate(jac->comm_hypre,&jac->hsolver); 79216d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut; 79316d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_Pilut; 79416d9e3a6SLisandro Dalcin jac->destroy = HYPRE_ParCSRPilutDestroy; 79516d9e3a6SLisandro Dalcin jac->setup = HYPRE_ParCSRPilutSetup; 79616d9e3a6SLisandro Dalcin jac->solve = HYPRE_ParCSRPilutSolve; 79716d9e3a6SLisandro Dalcin jac->factorrowsize = PETSC_DEFAULT; 79816d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 79916d9e3a6SLisandro Dalcin } 80016d9e3a6SLisandro Dalcin ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr); 80116d9e3a6SLisandro Dalcin if (flag) { 80216d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsCreate(jac->comm_hypre,&jac->hsolver); 80316d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails; 80416d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_ParaSails; 80516d9e3a6SLisandro Dalcin jac->destroy = HYPRE_ParaSailsDestroy; 80616d9e3a6SLisandro Dalcin jac->setup = HYPRE_ParaSailsSetup; 80716d9e3a6SLisandro Dalcin jac->solve = HYPRE_ParaSailsSolve; 80816d9e3a6SLisandro Dalcin /* initialize */ 80916d9e3a6SLisandro Dalcin jac->nlevels = 1; 81016d9e3a6SLisandro Dalcin jac->threshhold = .1; 81116d9e3a6SLisandro Dalcin jac->filter = .1; 81216d9e3a6SLisandro Dalcin jac->loadbal = 0; 81316d9e3a6SLisandro Dalcin if (PetscLogPrintInfo) { 81416d9e3a6SLisandro Dalcin jac->logging = (int) PETSC_TRUE; 81516d9e3a6SLisandro Dalcin } else { 81616d9e3a6SLisandro Dalcin jac->logging = (int) PETSC_FALSE; 81716d9e3a6SLisandro Dalcin } 81816d9e3a6SLisandro Dalcin jac->ruse = (int) PETSC_FALSE; 81916d9e3a6SLisandro Dalcin jac->symt = 0; 82016d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr); 82116d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr); 82216d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr); 82316d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr); 82416d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr); 82516d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr); 82616d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 82716d9e3a6SLisandro Dalcin } 82816d9e3a6SLisandro Dalcin ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr); 82916d9e3a6SLisandro Dalcin if (flag) { 83016d9e3a6SLisandro Dalcin ierr = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver); 83116d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid; 83216d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_Euclid; 83316d9e3a6SLisandro Dalcin jac->destroy = HYPRE_EuclidDestroy; 83416d9e3a6SLisandro Dalcin jac->setup = HYPRE_EuclidSetup; 83516d9e3a6SLisandro Dalcin jac->solve = HYPRE_EuclidSolve; 83616d9e3a6SLisandro Dalcin /* initialization */ 83716d9e3a6SLisandro Dalcin jac->bjilu = PETSC_FALSE; 83816d9e3a6SLisandro Dalcin jac->levels = 1; 83916d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 84016d9e3a6SLisandro Dalcin } 84116d9e3a6SLisandro Dalcin ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr); 84216d9e3a6SLisandro Dalcin if (flag) { 84316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGCreate(&jac->hsolver); 84416d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_BoomerAMG; 84516d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_BoomerAMG; 84616d9e3a6SLisandro Dalcin pc->ops->applytranspose = PCApplyTranspose_HYPRE_BoomerAMG; 84716d9e3a6SLisandro Dalcin pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG; 84816d9e3a6SLisandro Dalcin jac->destroy = HYPRE_BoomerAMGDestroy; 84916d9e3a6SLisandro Dalcin jac->setup = HYPRE_BoomerAMGSetup; 85016d9e3a6SLisandro Dalcin jac->solve = HYPRE_BoomerAMGSolve; 85116d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_FALSE; 85216d9e3a6SLisandro Dalcin /* these defaults match the hypre defaults */ 85316d9e3a6SLisandro Dalcin jac->cycletype = 1; 85416d9e3a6SLisandro Dalcin jac->maxlevels = 25; 85516d9e3a6SLisandro Dalcin jac->maxiter = 1; 8568f87f92bSBarry Smith jac->tol = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */ 85716d9e3a6SLisandro Dalcin jac->truncfactor = 0.0; 85816d9e3a6SLisandro Dalcin jac->strongthreshold = .25; 85916d9e3a6SLisandro Dalcin jac->maxrowsum = .9; 86016d9e3a6SLisandro Dalcin jac->coarsentype = 6; 86116d9e3a6SLisandro Dalcin jac->measuretype = 0; 8620f1074feSSatish Balay jac->gridsweeps[0] = jac->gridsweeps[1] = jac->gridsweeps[2] = 1; 8638f87f92bSBarry Smith jac->relaxtype[0] = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */ 8640f1074feSSatish Balay jac->relaxtype[2] = 9; /*G.E. */ 86516d9e3a6SLisandro Dalcin jac->relaxweight = 1.0; 86616d9e3a6SLisandro Dalcin jac->outerrelaxweight = 1.0; 86716d9e3a6SLisandro Dalcin jac->relaxorder = 1; 8680f1074feSSatish Balay jac->interptype = 0; 8690f1074feSSatish Balay jac->agg_nl = 0; 8700f1074feSSatish Balay jac->pmax = 0; 8710f1074feSSatish Balay jac->truncfactor = 0.0; 8720f1074feSSatish Balay jac->agg_num_paths = 1; 8738f87f92bSBarry Smith 8748f87f92bSBarry Smith jac->nodal_coarsen = 0; 8758f87f92bSBarry Smith jac->nodal_relax = PETSC_FALSE; 8768f87f92bSBarry Smith jac->nodal_relax_levels = 1; 87716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr); 87816d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr); 87916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 88016d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 88116d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr); 88216d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr); 88316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr); 88416d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr); 88516d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr); 88616d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr); 8870f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);CHKERRQ(ierr); 8880f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl); 8890f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);CHKERRQ(ierr); 8900f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);CHKERRQ(ierr); 8910f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, jac->relaxtype[0]); /*defaults coarse to 9*/ 8920f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver, jac->gridsweeps[0]); /*defaults coarse to 1 */ 89316d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 89416d9e3a6SLisandro Dalcin } 89516d9e3a6SLisandro Dalcin ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr); 89616d9e3a6SLisandro Dalcin jac->hypre_type = PETSC_NULL; 89716d9e3a6SLisandro Dalcin SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name); 89816d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 89916d9e3a6SLisandro Dalcin } 90016d9e3a6SLisandro Dalcin EXTERN_C_END 90116d9e3a6SLisandro Dalcin 90216d9e3a6SLisandro Dalcin /* 90316d9e3a6SLisandro Dalcin It only gets here if the HYPRE type has not been set before the call to 90416d9e3a6SLisandro Dalcin ...SetFromOptions() which actually is most of the time 90516d9e3a6SLisandro Dalcin */ 90616d9e3a6SLisandro Dalcin #undef __FUNCT__ 90716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE" 90816d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE(PC pc) 90916d9e3a6SLisandro Dalcin { 91016d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 91116d9e3a6SLisandro Dalcin PetscErrorCode ierr; 91216d9e3a6SLisandro Dalcin int indx; 91316d9e3a6SLisandro Dalcin const char *type[] = {"pilut","parasails","boomeramg","euclid"}; 91416d9e3a6SLisandro Dalcin PetscTruth flg; 91516d9e3a6SLisandro Dalcin 91616d9e3a6SLisandro Dalcin PetscFunctionBegin; 91716d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr); 91816d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"pilut",&indx,&flg);CHKERRQ(ierr); 91916d9e3a6SLisandro Dalcin if (PetscOptionsPublishCount) { /* force the default if it was not yet set and user did not set with option */ 92016d9e3a6SLisandro Dalcin if (!flg && !jac->hypre_type) { 92116d9e3a6SLisandro Dalcin flg = PETSC_TRUE; 92216d9e3a6SLisandro Dalcin indx = 0; 92316d9e3a6SLisandro Dalcin } 92416d9e3a6SLisandro Dalcin } 92516d9e3a6SLisandro Dalcin if (flg) { 92616d9e3a6SLisandro Dalcin ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr); 92716d9e3a6SLisandro Dalcin } 92816d9e3a6SLisandro Dalcin if (pc->ops->setfromoptions) { 92916d9e3a6SLisandro Dalcin ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr); 93016d9e3a6SLisandro Dalcin } 93116d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 93216d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 93316d9e3a6SLisandro Dalcin } 93416d9e3a6SLisandro Dalcin 93516d9e3a6SLisandro Dalcin #undef __FUNCT__ 93616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType" 93716d9e3a6SLisandro Dalcin /*@C 93816d9e3a6SLisandro Dalcin PCHYPRESetType - Sets which hypre preconditioner you wish to use 93916d9e3a6SLisandro Dalcin 94016d9e3a6SLisandro Dalcin Input Parameters: 94116d9e3a6SLisandro Dalcin + pc - the preconditioner context 94216d9e3a6SLisandro Dalcin - name - either pilut, parasails, boomeramg, euclid 94316d9e3a6SLisandro Dalcin 94416d9e3a6SLisandro Dalcin Options Database Keys: 94516d9e3a6SLisandro Dalcin -pc_hypre_type - One of pilut, parasails, boomeramg, euclid 94616d9e3a6SLisandro Dalcin 94716d9e3a6SLisandro Dalcin Level: intermediate 94816d9e3a6SLisandro Dalcin 94916d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 95016d9e3a6SLisandro Dalcin PCHYPRE 95116d9e3a6SLisandro Dalcin 95216d9e3a6SLisandro Dalcin @*/ 95316d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC pc,const char name[]) 95416d9e3a6SLisandro Dalcin { 95516d9e3a6SLisandro Dalcin PetscErrorCode ierr,(*f)(PC,const char[]); 95616d9e3a6SLisandro Dalcin 95716d9e3a6SLisandro Dalcin PetscFunctionBegin; 95816d9e3a6SLisandro Dalcin PetscValidHeaderSpecific(pc,PC_COOKIE,1); 95916d9e3a6SLisandro Dalcin PetscValidCharPointer(name,2); 96016d9e3a6SLisandro Dalcin ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPRESetType_C",(void (**)(void))&f);CHKERRQ(ierr); 96116d9e3a6SLisandro Dalcin if (f) { 96216d9e3a6SLisandro Dalcin ierr = (*f)(pc,name);CHKERRQ(ierr); 96316d9e3a6SLisandro Dalcin } 96416d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 96516d9e3a6SLisandro Dalcin } 96616d9e3a6SLisandro Dalcin 96716d9e3a6SLisandro Dalcin #undef __FUNCT__ 96816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType" 96916d9e3a6SLisandro Dalcin /*@C 97016d9e3a6SLisandro Dalcin PCHYPREGetType - Gets which hypre preconditioner you are using 97116d9e3a6SLisandro Dalcin 97216d9e3a6SLisandro Dalcin Input Parameter: 97316d9e3a6SLisandro Dalcin . pc - the preconditioner context 97416d9e3a6SLisandro Dalcin 97516d9e3a6SLisandro Dalcin Output Parameter: 97616d9e3a6SLisandro Dalcin . name - either pilut, parasails, boomeramg, euclid 97716d9e3a6SLisandro Dalcin 97816d9e3a6SLisandro Dalcin Level: intermediate 97916d9e3a6SLisandro Dalcin 98016d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC, 98116d9e3a6SLisandro Dalcin PCHYPRE 98216d9e3a6SLisandro Dalcin 98316d9e3a6SLisandro Dalcin @*/ 98416d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC pc,const char *name[]) 98516d9e3a6SLisandro Dalcin { 98616d9e3a6SLisandro Dalcin PetscErrorCode ierr,(*f)(PC,const char*[]); 98716d9e3a6SLisandro Dalcin 98816d9e3a6SLisandro Dalcin PetscFunctionBegin; 98916d9e3a6SLisandro Dalcin PetscValidHeaderSpecific(pc,PC_COOKIE,1); 99016d9e3a6SLisandro Dalcin PetscValidPointer(name,2); 99116d9e3a6SLisandro Dalcin ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPREGetType_C",(void (**)(void))&f);CHKERRQ(ierr); 99216d9e3a6SLisandro Dalcin if (f) { 99316d9e3a6SLisandro Dalcin ierr = (*f)(pc,name);CHKERRQ(ierr); 99416d9e3a6SLisandro Dalcin } 99516d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 99616d9e3a6SLisandro Dalcin } 99716d9e3a6SLisandro Dalcin 99816d9e3a6SLisandro Dalcin /*MC 99916d9e3a6SLisandro Dalcin PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre 100016d9e3a6SLisandro Dalcin 100116d9e3a6SLisandro Dalcin Options Database Keys: 100216d9e3a6SLisandro Dalcin + -pc_hypre_type - One of pilut, parasails, boomeramg, euclid 100316d9e3a6SLisandro Dalcin - Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX 100416d9e3a6SLisandro Dalcin preconditioner 100516d9e3a6SLisandro Dalcin 100616d9e3a6SLisandro Dalcin Level: intermediate 100716d9e3a6SLisandro Dalcin 100816d9e3a6SLisandro Dalcin Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()), 100916d9e3a6SLisandro Dalcin the many hypre options can ONLY be set via the options database (e.g. the command line 101016d9e3a6SLisandro Dalcin or with PetscOptionsSetValue(), there are no functions to set them) 101116d9e3a6SLisandro Dalcin 101216d9e3a6SLisandro Dalcin The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations 10130f1074feSSatish Balay (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if 10140f1074feSSatish Balay -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner 10150f1074feSSatish Balay (-pc_hypre_boomeramg_rtol should be set to 0.0 - the default - to strictly use a fixed number of 10168f87f92bSBarry Smith iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations 10170f1074feSSatish Balay and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10 10180f1074feSSatish Balay then AT MOST twenty V-cycles of boomeramg will be called. 101916d9e3a6SLisandro Dalcin 10200f1074feSSatish Balay Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation 10210f1074feSSatish Balay (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry. 10220f1074feSSatish Balay Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi. 102316d9e3a6SLisandro Dalcin If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly 102416d9e3a6SLisandro Dalcin and use -ksp_max_it to control the number of V-cycles. 102516d9e3a6SLisandro Dalcin (see the PETSc FAQ.html at the PETSc website under the Documentation tab). 102616d9e3a6SLisandro Dalcin 102716d9e3a6SLisandro Dalcin 2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option 102816d9e3a6SLisandro Dalcin -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L. 102916d9e3a6SLisandro Dalcin 103016d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 103116d9e3a6SLisandro Dalcin PCHYPRESetType() 103216d9e3a6SLisandro Dalcin 103316d9e3a6SLisandro Dalcin M*/ 103416d9e3a6SLisandro Dalcin 103516d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 103616d9e3a6SLisandro Dalcin #undef __FUNCT__ 103716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCCreate_HYPRE" 103816d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_HYPRE(PC pc) 103916d9e3a6SLisandro Dalcin { 104016d9e3a6SLisandro Dalcin PC_HYPRE *jac; 104116d9e3a6SLisandro Dalcin PetscErrorCode ierr; 104216d9e3a6SLisandro Dalcin 104316d9e3a6SLisandro Dalcin PetscFunctionBegin; 1044*38f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_HYPRE,&jac);CHKERRQ(ierr); 104516d9e3a6SLisandro Dalcin pc->data = jac; 104616d9e3a6SLisandro Dalcin pc->ops->destroy = PCDestroy_HYPRE; 104716d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE; 104816d9e3a6SLisandro Dalcin pc->ops->setup = PCSetUp_HYPRE; 104916d9e3a6SLisandro Dalcin pc->ops->apply = PCApply_HYPRE; 105016d9e3a6SLisandro Dalcin jac->comm_hypre = MPI_COMM_NULL; 105116d9e3a6SLisandro Dalcin jac->hypre_type = PETSC_NULL; 105216d9e3a6SLisandro Dalcin /* duplicate communicator for hypre */ 105316d9e3a6SLisandro Dalcin ierr = MPI_Comm_dup(pc->comm,&(jac->comm_hypre));CHKERRQ(ierr); 105416d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr); 105516d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr); 105616d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 105716d9e3a6SLisandro Dalcin } 105816d9e3a6SLisandro Dalcin EXTERN_C_END 1059