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); 196745373dSBarry Smith EXTERN PetscErrorCode MatHYPRE_IJMatrixFastCopy(Mat,HYPRE_IJMatrix); 2016d9e3a6SLisandro Dalcin EXTERN PetscErrorCode VecHYPRE_IJVectorCreate(Vec,HYPRE_IJVector*); 2116d9e3a6SLisandro Dalcin 2216d9e3a6SLisandro Dalcin /* 2316d9e3a6SLisandro Dalcin Private context (data structure) for the preconditioner. 2416d9e3a6SLisandro Dalcin */ 2516d9e3a6SLisandro Dalcin typedef struct { 2616d9e3a6SLisandro Dalcin HYPRE_Solver hsolver; 2716d9e3a6SLisandro Dalcin HYPRE_IJMatrix ij; 2816d9e3a6SLisandro Dalcin HYPRE_IJVector b,x; 2916d9e3a6SLisandro Dalcin 3016d9e3a6SLisandro Dalcin PetscErrorCode (*destroy)(HYPRE_Solver); 3116d9e3a6SLisandro Dalcin PetscErrorCode (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector); 3216d9e3a6SLisandro Dalcin PetscErrorCode (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector); 3316d9e3a6SLisandro Dalcin 3416d9e3a6SLisandro Dalcin MPI_Comm comm_hypre; 3516d9e3a6SLisandro Dalcin char *hypre_type; 3616d9e3a6SLisandro Dalcin 3716d9e3a6SLisandro Dalcin /* options for Pilut and BoomerAMG*/ 3816d9e3a6SLisandro Dalcin int maxiter; 3916d9e3a6SLisandro Dalcin double tol; 4016d9e3a6SLisandro Dalcin 4116d9e3a6SLisandro Dalcin /* options for Pilut */ 4216d9e3a6SLisandro Dalcin int factorrowsize; 4316d9e3a6SLisandro Dalcin 4416d9e3a6SLisandro Dalcin /* options for ParaSails */ 4516d9e3a6SLisandro Dalcin int nlevels; 4616d9e3a6SLisandro Dalcin double threshhold; 4716d9e3a6SLisandro Dalcin double filter; 4816d9e3a6SLisandro Dalcin int sym; 4916d9e3a6SLisandro Dalcin double loadbal; 5016d9e3a6SLisandro Dalcin int logging; 5116d9e3a6SLisandro Dalcin int ruse; 5216d9e3a6SLisandro Dalcin int symt; 5316d9e3a6SLisandro Dalcin 5416d9e3a6SLisandro Dalcin /* options for Euclid */ 5516d9e3a6SLisandro Dalcin PetscTruth bjilu; 5616d9e3a6SLisandro Dalcin int levels; 5716d9e3a6SLisandro Dalcin 5816d9e3a6SLisandro Dalcin /* options for Euclid and BoomerAMG */ 5916d9e3a6SLisandro Dalcin PetscTruth printstatistics; 6016d9e3a6SLisandro Dalcin 6116d9e3a6SLisandro Dalcin /* options for BoomerAMG */ 6216d9e3a6SLisandro Dalcin int cycletype; 6316d9e3a6SLisandro Dalcin int maxlevels; 6416d9e3a6SLisandro Dalcin double strongthreshold; 6516d9e3a6SLisandro Dalcin double maxrowsum; 660f1074feSSatish Balay int gridsweeps[3]; 6716d9e3a6SLisandro Dalcin int coarsentype; 6816d9e3a6SLisandro Dalcin int measuretype; 690f1074feSSatish Balay int relaxtype[3]; 7016d9e3a6SLisandro Dalcin double relaxweight; 7116d9e3a6SLisandro Dalcin double outerrelaxweight; 7216d9e3a6SLisandro Dalcin int relaxorder; 7316d9e3a6SLisandro Dalcin double truncfactor; 7416d9e3a6SLisandro Dalcin PetscTruth applyrichardson; 750f1074feSSatish Balay int pmax; 760f1074feSSatish Balay int interptype; 770f1074feSSatish Balay int agg_nl; 780f1074feSSatish Balay int agg_num_paths; 798f87f92bSBarry Smith int nodal_coarsen; 808f87f92bSBarry Smith PetscTruth nodal_relax; 818f87f92bSBarry Smith int nodal_relax_levels; 8216d9e3a6SLisandro Dalcin } PC_HYPRE; 8316d9e3a6SLisandro Dalcin 8416d9e3a6SLisandro Dalcin 8516d9e3a6SLisandro Dalcin #undef __FUNCT__ 8616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetUp_HYPRE" 8716d9e3a6SLisandro Dalcin static PetscErrorCode PCSetUp_HYPRE(PC pc) 8816d9e3a6SLisandro Dalcin { 8916d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 9016d9e3a6SLisandro Dalcin PetscErrorCode ierr; 9116d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 9216d9e3a6SLisandro Dalcin HYPRE_ParVector bv,xv; 9316d9e3a6SLisandro Dalcin PetscInt bs; 9416d9e3a6SLisandro Dalcin int hierr; 9516d9e3a6SLisandro Dalcin 9616d9e3a6SLisandro Dalcin PetscFunctionBegin; 9716d9e3a6SLisandro Dalcin if (!jac->hypre_type) { 9802a17cd4SBarry Smith ierr = PCHYPRESetType(pc,"boomeramg");CHKERRQ(ierr); 9916d9e3a6SLisandro Dalcin } 1005f5c5b43SBarry Smith 1015f5c5b43SBarry Smith if (pc->setupcalled) { 1025f5c5b43SBarry Smith /* always destroy the old matrix and create a new memory; 1035f5c5b43SBarry Smith hope this does not churn the memory too much. The problem 1045f5c5b43SBarry Smith is I do not know if it is possible to put the matrix back to 1055f5c5b43SBarry Smith its initial state so that we can directly copy the values 1065f5c5b43SBarry Smith the second time through. */ 10716d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); 1085f5c5b43SBarry Smith jac->ij = 0; 10916d9e3a6SLisandro Dalcin } 1105f5c5b43SBarry Smith 11116d9e3a6SLisandro Dalcin if (!jac->ij) { /* create the matrix the first time through */ 11216d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr); 11316d9e3a6SLisandro Dalcin } 11416d9e3a6SLisandro Dalcin if (!jac->b) { /* create the vectors the first time through */ 11516d9e3a6SLisandro Dalcin Vec x,b; 11616d9e3a6SLisandro Dalcin ierr = MatGetVecs(pc->pmat,&x,&b);CHKERRQ(ierr); 11716d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr); 11816d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr); 11916d9e3a6SLisandro Dalcin ierr = VecDestroy(x);CHKERRQ(ierr); 12016d9e3a6SLisandro Dalcin ierr = VecDestroy(b);CHKERRQ(ierr); 12116d9e3a6SLisandro Dalcin } 1225f5c5b43SBarry Smith 12316d9e3a6SLisandro Dalcin /* special case for BoomerAMG */ 12416d9e3a6SLisandro Dalcin if (jac->setup == HYPRE_BoomerAMGSetup) { 12516d9e3a6SLisandro Dalcin ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr); 12616d9e3a6SLisandro Dalcin if (bs > 1) { 12716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetNumFunctions(jac->hsolver,bs);CHKERRQ(ierr); 12816d9e3a6SLisandro Dalcin } 12916d9e3a6SLisandro Dalcin }; 13016d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr); 13116d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 13216d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&bv);CHKERRQ(ierr); 13316d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&xv);CHKERRQ(ierr); 13416d9e3a6SLisandro Dalcin hierr = (*jac->setup)(jac->hsolver,hmat,bv,xv); 13516d9e3a6SLisandro Dalcin if (hierr) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE setup, error code %d",hierr); 13616d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 13716d9e3a6SLisandro Dalcin } 13816d9e3a6SLisandro Dalcin 13916d9e3a6SLisandro Dalcin /* 14016d9e3a6SLisandro Dalcin Replaces the address where the HYPRE vector points to its data with the address of 14116d9e3a6SLisandro Dalcin PETSc's data. Saves the old address so it can be reset when we are finished with it. 14216d9e3a6SLisandro Dalcin Allows use to get the data into a HYPRE vector without the cost of memcopies 14316d9e3a6SLisandro Dalcin */ 14416d9e3a6SLisandro Dalcin #define HYPREReplacePointer(b,newvalue,savedvalue) {\ 14516d9e3a6SLisandro Dalcin hypre_ParVector *par_vector = (hypre_ParVector *)hypre_IJVectorObject(((hypre_IJVector*)b));\ 14616d9e3a6SLisandro Dalcin hypre_Vector *local_vector = hypre_ParVectorLocalVector(par_vector);\ 14716d9e3a6SLisandro Dalcin savedvalue = local_vector->data;\ 14816d9e3a6SLisandro Dalcin local_vector->data = newvalue;} 14916d9e3a6SLisandro Dalcin 15016d9e3a6SLisandro Dalcin #undef __FUNCT__ 15116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApply_HYPRE" 15216d9e3a6SLisandro Dalcin static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x) 15316d9e3a6SLisandro Dalcin { 15416d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 15516d9e3a6SLisandro Dalcin PetscErrorCode ierr; 15616d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 15716d9e3a6SLisandro Dalcin PetscScalar *bv,*xv; 15816d9e3a6SLisandro Dalcin HYPRE_ParVector jbv,jxv; 15916d9e3a6SLisandro Dalcin PetscScalar *sbv,*sxv; 16016d9e3a6SLisandro Dalcin int hierr; 16116d9e3a6SLisandro Dalcin 16216d9e3a6SLisandro Dalcin PetscFunctionBegin; 16316d9e3a6SLisandro Dalcin if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);} 16416d9e3a6SLisandro Dalcin ierr = VecGetArray(b,&bv);CHKERRQ(ierr); 16516d9e3a6SLisandro Dalcin ierr = VecGetArray(x,&xv);CHKERRQ(ierr); 16616d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,bv,sbv); 16716d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,xv,sxv); 16816d9e3a6SLisandro Dalcin 16916d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 17016d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr); 17116d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr); 17216d9e3a6SLisandro Dalcin hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv); 1730f1074feSSatish Balay 1740f1074feSSatish Balay /*if (hierr && (hierr != HYPRE_ERROR_CONV || jac->solve != HYPRE_BoomerAMGSolve))SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 1750f1074feSSatish Balay */ 1760f1074feSSatish Balay /* error code of HYPRE_ERROR_CONV means convergence not achieved - if 1770f1074feSSatish Balay the tolerance is set to 0.0 (the default), a convergence error will 1780f1074feSSatish Balay not occur (so we may not want to overide the conv. error here?*/ 1790f1074feSSatish Balay if (hierr && hierr != HYPRE_ERROR_CONV) 1800f1074feSSatish Balay { 1810f1074feSSatish Balay SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 1820f1074feSSatish Balay } 18316d9e3a6SLisandro Dalcin if (hierr) hypre__global_error = 0; 18416d9e3a6SLisandro Dalcin 1850f1074feSSatish Balay 18616d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,sbv,bv); 18716d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,sxv,xv); 18816d9e3a6SLisandro Dalcin ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr); 18916d9e3a6SLisandro Dalcin ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr); 19016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 19116d9e3a6SLisandro Dalcin } 19216d9e3a6SLisandro Dalcin 19316d9e3a6SLisandro Dalcin #undef __FUNCT__ 19416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCDestroy_HYPRE" 19516d9e3a6SLisandro Dalcin static PetscErrorCode PCDestroy_HYPRE(PC pc) 19616d9e3a6SLisandro Dalcin { 19716d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 19816d9e3a6SLisandro Dalcin PetscErrorCode ierr; 19916d9e3a6SLisandro Dalcin 20016d9e3a6SLisandro Dalcin PetscFunctionBegin; 20116d9e3a6SLisandro Dalcin if (jac->ij) { ierr = HYPRE_IJMatrixDestroy(jac->ij);CHKERRQ(ierr); } 20216d9e3a6SLisandro Dalcin if (jac->b) { ierr = HYPRE_IJVectorDestroy(jac->b);CHKERRQ(ierr); } 20316d9e3a6SLisandro Dalcin if (jac->x) { ierr = HYPRE_IJVectorDestroy(jac->x);CHKERRQ(ierr); } 20416d9e3a6SLisandro Dalcin if (jac->destroy) { ierr = (*jac->destroy)(jac->hsolver);CHKERRQ(ierr); } 20516d9e3a6SLisandro Dalcin ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr); 20616d9e3a6SLisandro Dalcin if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);} 20716d9e3a6SLisandro Dalcin ierr = PetscFree(jac);CHKERRQ(ierr); 20816d9e3a6SLisandro Dalcin 20916d9e3a6SLisandro Dalcin ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr); 21016d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","",PETSC_NULL);CHKERRQ(ierr); 21116d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","",PETSC_NULL);CHKERRQ(ierr); 21216d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 21316d9e3a6SLisandro Dalcin } 21416d9e3a6SLisandro Dalcin 21516d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 21616d9e3a6SLisandro Dalcin #undef __FUNCT__ 21716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut" 21816d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PC pc) 21916d9e3a6SLisandro Dalcin { 22016d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 22116d9e3a6SLisandro Dalcin PetscErrorCode ierr; 22216d9e3a6SLisandro Dalcin PetscTruth flag; 22316d9e3a6SLisandro Dalcin 22416d9e3a6SLisandro Dalcin PetscFunctionBegin; 22516d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE Pilut Options");CHKERRQ(ierr); 22616d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr); 22716d9e3a6SLisandro Dalcin if (flag) { 22816d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 22916d9e3a6SLisandro Dalcin } 23016d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr); 23116d9e3a6SLisandro Dalcin if (flag) { 23216d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutSetDropTolerance(jac->hsolver,jac->tol);CHKERRQ(ierr); 23316d9e3a6SLisandro Dalcin } 23416d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr); 23516d9e3a6SLisandro Dalcin if (flag) { 23616d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutSetFactorRowSize(jac->hsolver,jac->factorrowsize);CHKERRQ(ierr); 23716d9e3a6SLisandro Dalcin } 23816d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 23916d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 24016d9e3a6SLisandro Dalcin } 24116d9e3a6SLisandro Dalcin 24216d9e3a6SLisandro Dalcin #undef __FUNCT__ 24316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Pilut" 24416d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer) 24516d9e3a6SLisandro Dalcin { 24616d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 24716d9e3a6SLisandro Dalcin PetscErrorCode ierr; 24816d9e3a6SLisandro Dalcin PetscTruth iascii; 24916d9e3a6SLisandro Dalcin 25016d9e3a6SLisandro Dalcin PetscFunctionBegin; 25116d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 25216d9e3a6SLisandro Dalcin if (iascii) { 25316d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut preconditioning\n");CHKERRQ(ierr); 25416d9e3a6SLisandro Dalcin if (jac->maxiter != PETSC_DEFAULT) { 25516d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr); 25616d9e3a6SLisandro Dalcin } else { 25716d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr); 25816d9e3a6SLisandro Dalcin } 25916d9e3a6SLisandro Dalcin if (jac->tol != PETSC_DEFAULT) { 26016d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: drop tolerance %G\n",jac->tol);CHKERRQ(ierr); 26116d9e3a6SLisandro Dalcin } else { 26216d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr); 26316d9e3a6SLisandro Dalcin } 26416d9e3a6SLisandro Dalcin if (jac->factorrowsize != PETSC_DEFAULT) { 26516d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr); 26616d9e3a6SLisandro Dalcin } else { 26716d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default factor row size \n");CHKERRQ(ierr); 26816d9e3a6SLisandro Dalcin } 26916d9e3a6SLisandro Dalcin } 27016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 27116d9e3a6SLisandro Dalcin } 27216d9e3a6SLisandro Dalcin 27316d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 27416d9e3a6SLisandro Dalcin #undef __FUNCT__ 27516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Euclid" 27616d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_Euclid(PC pc) 27716d9e3a6SLisandro Dalcin { 27816d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 27916d9e3a6SLisandro Dalcin PetscErrorCode ierr; 28016d9e3a6SLisandro Dalcin PetscTruth flag; 281390e7148SBarry Smith char *args[8],levels[16]; 282390e7148SBarry Smith PetscInt cnt = 0; 28316d9e3a6SLisandro Dalcin 28416d9e3a6SLisandro Dalcin PetscFunctionBegin; 28516d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE Euclid Options");CHKERRQ(ierr); 28616d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_euclid_levels","Number of levels of fill ILU(k)","None",jac->levels,&jac->levels,&flag);CHKERRQ(ierr); 28716d9e3a6SLisandro Dalcin if (flag) { 28816d9e3a6SLisandro Dalcin if (jac->levels < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be nonegative",jac->levels); 28916d9e3a6SLisandro Dalcin sprintf(levels,"%d",jac->levels); 290390e7148SBarry Smith args[cnt++] = (char*)"-level"; args[cnt++] = levels; 29116d9e3a6SLisandro Dalcin } 29216d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_euclid_bj","Use block Jacobi ILU(k)","None",jac->bjilu,&jac->bjilu,PETSC_NULL);CHKERRQ(ierr); 29316d9e3a6SLisandro Dalcin if (jac->bjilu) { 294390e7148SBarry Smith args[cnt++] =(char*) "-bj"; args[cnt++] = (char*)"1"; 29516d9e3a6SLisandro Dalcin } 29616d9e3a6SLisandro Dalcin 29716d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_euclid_print_statistics","Print statistics","None",jac->printstatistics,&jac->printstatistics,PETSC_NULL);CHKERRQ(ierr); 29816d9e3a6SLisandro Dalcin if (jac->printstatistics) { 299390e7148SBarry Smith args[cnt++] = (char*)"-eu_stats"; args[cnt++] = (char*)"1"; 300390e7148SBarry Smith args[cnt++] = (char*)"-eu_mem"; args[cnt++] = (char*)"1"; 30116d9e3a6SLisandro Dalcin } 30216d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 303390e7148SBarry Smith if (cnt) { 304390e7148SBarry Smith ierr = HYPRE_EuclidSetParams(jac->hsolver,cnt,args);CHKERRQ(ierr); 305390e7148SBarry Smith } 30616d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 30716d9e3a6SLisandro Dalcin } 30816d9e3a6SLisandro Dalcin 30916d9e3a6SLisandro Dalcin #undef __FUNCT__ 31016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Euclid" 31116d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Euclid(PC pc,PetscViewer viewer) 31216d9e3a6SLisandro Dalcin { 31316d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 31416d9e3a6SLisandro Dalcin PetscErrorCode ierr; 31516d9e3a6SLisandro Dalcin PetscTruth iascii; 31616d9e3a6SLisandro Dalcin 31716d9e3a6SLisandro Dalcin PetscFunctionBegin; 31816d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 31916d9e3a6SLisandro Dalcin if (iascii) { 32016d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid preconditioning\n");CHKERRQ(ierr); 32116d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid: number of levels %d\n",jac->levels);CHKERRQ(ierr); 32216d9e3a6SLisandro Dalcin if (jac->bjilu) { 32316d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Euclid: Using block Jacobi ILU instead of parallel ILU\n");CHKERRQ(ierr); 32416d9e3a6SLisandro Dalcin } 32516d9e3a6SLisandro Dalcin } 32616d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 32716d9e3a6SLisandro Dalcin } 32816d9e3a6SLisandro Dalcin 32916d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 33016d9e3a6SLisandro Dalcin 33116d9e3a6SLisandro Dalcin #undef __FUNCT__ 33216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG" 33316d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x) 33416d9e3a6SLisandro Dalcin { 33516d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 33616d9e3a6SLisandro Dalcin PetscErrorCode ierr; 33716d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 33816d9e3a6SLisandro Dalcin PetscScalar *bv,*xv; 33916d9e3a6SLisandro Dalcin HYPRE_ParVector jbv,jxv; 34016d9e3a6SLisandro Dalcin PetscScalar *sbv,*sxv; 34116d9e3a6SLisandro Dalcin int hierr; 34216d9e3a6SLisandro Dalcin 34316d9e3a6SLisandro Dalcin PetscFunctionBegin; 34416d9e3a6SLisandro Dalcin ierr = VecSet(x,0.0);CHKERRQ(ierr); 34516d9e3a6SLisandro Dalcin ierr = VecGetArray(b,&bv);CHKERRQ(ierr); 34616d9e3a6SLisandro Dalcin ierr = VecGetArray(x,&xv);CHKERRQ(ierr); 34716d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,bv,sbv); 34816d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,xv,sxv); 34916d9e3a6SLisandro Dalcin 35016d9e3a6SLisandro Dalcin ierr = HYPRE_IJMatrixGetObject(jac->ij,(void**)&hmat);CHKERRQ(ierr); 35116d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->b,(void**)&jbv);CHKERRQ(ierr); 35216d9e3a6SLisandro Dalcin ierr = HYPRE_IJVectorGetObject(jac->x,(void**)&jxv);CHKERRQ(ierr); 35316d9e3a6SLisandro Dalcin 35416d9e3a6SLisandro Dalcin hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv); 35516d9e3a6SLisandro Dalcin /* error code of 1 in BoomerAMG merely means convergence not achieved */ 35616d9e3a6SLisandro Dalcin if (hierr && (hierr != 1)) SETERRQ1(PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 35716d9e3a6SLisandro Dalcin if (hierr) hypre__global_error = 0; 35816d9e3a6SLisandro Dalcin 35916d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,sbv,bv); 36016d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,sxv,xv); 36116d9e3a6SLisandro Dalcin ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr); 36216d9e3a6SLisandro Dalcin ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr); 36316d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 36416d9e3a6SLisandro Dalcin } 36516d9e3a6SLisandro Dalcin 36616d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCycleType[] = {"","V","W"}; 3670f1074feSSatish Balay static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"}; 36816d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGMeasureType[] = {"local","global"}; 3690f1074feSSatish Balay static const char *HYPREBoomerAMGRelaxType[] = {"Jacobi","sequential-Gauss-Seidel","","SOR/Jacobi","backward-SOR/Jacobi","","symmetric-SOR/Jacobi", 3700f1074feSSatish Balay "","","Gaussian-elimination"}; 3710f1074feSSatish Balay static const char *HYPREBoomerAMGInterpType[] = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i", 3720f1074feSSatish Balay "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"}; 37316d9e3a6SLisandro Dalcin #undef __FUNCT__ 37416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG" 37516d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PC pc) 37616d9e3a6SLisandro Dalcin { 37716d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 37816d9e3a6SLisandro Dalcin PetscErrorCode ierr; 37916d9e3a6SLisandro Dalcin int n,indx; 38016d9e3a6SLisandro Dalcin PetscTruth flg, tmp_truth; 38116d9e3a6SLisandro Dalcin double tmpdbl, twodbl[2]; 38216d9e3a6SLisandro Dalcin 38316d9e3a6SLisandro Dalcin PetscFunctionBegin; 38416d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE BoomerAMG Options");CHKERRQ(ierr); 38516d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr); 38616d9e3a6SLisandro Dalcin if (flg) { 38716d9e3a6SLisandro Dalcin jac->cycletype = indx; 38816d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr); 38916d9e3a6SLisandro Dalcin } 39016d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr); 39116d9e3a6SLisandro Dalcin if (flg) { 39216d9e3a6SLisandro Dalcin if (jac->maxlevels < 2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels); 39316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr); 39416d9e3a6SLisandro Dalcin } 39516d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr); 39616d9e3a6SLisandro Dalcin if (flg) { 39716d9e3a6SLisandro Dalcin if (jac->maxiter < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter); 39816d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 39916d9e3a6SLisandro Dalcin } 4000f1074feSSatish 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); 40116d9e3a6SLisandro Dalcin if (flg) { 4020f1074feSSatish Balay if (jac->tol < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %G must be greater than or equal to zero",jac->tol); 40316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 40416d9e3a6SLisandro Dalcin } 40516d9e3a6SLisandro Dalcin 4060f1074feSSatish Balay ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr); 40716d9e3a6SLisandro Dalcin if (flg) { 40816d9e3a6SLisandro Dalcin if (jac->truncfactor < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %G must be great than or equal zero",jac->truncfactor); 40916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr); 41016d9e3a6SLisandro Dalcin } 41116d9e3a6SLisandro Dalcin 4120f1074feSSatish 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); 4130f1074feSSatish Balay if (flg) { 4140f1074feSSatish Balay if (jac->pmax < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"P_max %G must be greater than or equal to zero",jac->pmax); 4150f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);CHKERRQ(ierr); 4160f1074feSSatish Balay } 4170f1074feSSatish Balay 4180f1074feSSatish Balay ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);CHKERRQ(ierr); 4190f1074feSSatish Balay if (flg) { 4200f1074feSSatish 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); 4210f1074feSSatish Balay 4220f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl);CHKERRQ(ierr); 4230f1074feSSatish Balay } 4240f1074feSSatish Balay 4250f1074feSSatish Balay 4260f1074feSSatish 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); 4270f1074feSSatish Balay if (flg) { 4280f1074feSSatish 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); 4290f1074feSSatish Balay 4300f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);CHKERRQ(ierr); 4310f1074feSSatish Balay } 4320f1074feSSatish Balay 4330f1074feSSatish Balay 43416d9e3a6SLisandro Dalcin ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr); 43516d9e3a6SLisandro Dalcin if (flg) { 43616d9e3a6SLisandro Dalcin if (jac->strongthreshold < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %G must be great than or equal zero",jac->strongthreshold); 43716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr); 43816d9e3a6SLisandro Dalcin } 43916d9e3a6SLisandro Dalcin ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr); 44016d9e3a6SLisandro Dalcin if (flg) { 44116d9e3a6SLisandro Dalcin if (jac->maxrowsum < 0.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be greater than zero",jac->maxrowsum); 44216d9e3a6SLisandro Dalcin if (jac->maxrowsum > 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %G must be less than or equal one",jac->maxrowsum); 44316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr); 44416d9e3a6SLisandro Dalcin } 44516d9e3a6SLisandro Dalcin 44616d9e3a6SLisandro Dalcin /* Grid sweeps */ 4470f1074feSSatish 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); 44816d9e3a6SLisandro Dalcin if (flg) { 44916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver,indx);CHKERRQ(ierr); 45016d9e3a6SLisandro Dalcin /* modify the jac structure so we can view the updated options with PC_View */ 45116d9e3a6SLisandro Dalcin jac->gridsweeps[0] = indx; 4520f1074feSSatish Balay jac->gridsweeps[1] = indx; 4530f1074feSSatish Balay /*defaults coarse to 1 */ 4540f1074feSSatish Balay jac->gridsweeps[2] = 1; 45516d9e3a6SLisandro Dalcin } 4560f1074feSSatish Balay 4570f1074feSSatish Balay ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx ,&flg);CHKERRQ(ierr); 45816d9e3a6SLisandro Dalcin if (flg) { 45916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 1);CHKERRQ(ierr); 4600f1074feSSatish Balay jac->gridsweeps[0] = indx; 46116d9e3a6SLisandro Dalcin } 46216d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr); 46316d9e3a6SLisandro Dalcin if (flg) { 46416d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 2);CHKERRQ(ierr); 4650f1074feSSatish Balay jac->gridsweeps[1] = indx; 46616d9e3a6SLisandro Dalcin } 4670f1074feSSatish Balay ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr); 46816d9e3a6SLisandro Dalcin if (flg) { 46916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleNumSweeps(jac->hsolver,indx, 3);CHKERRQ(ierr); 4700f1074feSSatish Balay jac->gridsweeps[2] = indx; 47116d9e3a6SLisandro Dalcin } 47216d9e3a6SLisandro Dalcin 47316d9e3a6SLisandro Dalcin /* Relax type */ 4740f1074feSSatish 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); 47516d9e3a6SLisandro Dalcin if (flg) { 4760f1074feSSatish Balay jac->relaxtype[0] = jac->relaxtype[1] = indx; 47716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, indx);CHKERRQ(ierr); 4780f1074feSSatish Balay /* by default, coarse type set to 9 */ 4790f1074feSSatish Balay jac->relaxtype[2] = 9; 4800f1074feSSatish Balay 48116d9e3a6SLisandro Dalcin } 4820f1074feSSatish Balay ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr); 48316d9e3a6SLisandro Dalcin if (flg) { 48416d9e3a6SLisandro Dalcin jac->relaxtype[0] = indx; 48516d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 1);CHKERRQ(ierr); 48616d9e3a6SLisandro Dalcin } 4870f1074feSSatish Balay ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr); 48816d9e3a6SLisandro Dalcin if (flg) { 4890f1074feSSatish Balay jac->relaxtype[1] = indx; 49016d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 2);CHKERRQ(ierr); 49116d9e3a6SLisandro Dalcin } 49216d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,10,HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr); 49316d9e3a6SLisandro Dalcin if (flg) { 4940f1074feSSatish Balay jac->relaxtype[2] = indx; 49516d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleRelaxType(jac->hsolver, indx, 3);CHKERRQ(ierr); 49616d9e3a6SLisandro Dalcin } 49716d9e3a6SLisandro Dalcin 49816d9e3a6SLisandro Dalcin /* Relaxation Weight */ 49916d9e3a6SLisandro 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); 50016d9e3a6SLisandro Dalcin if (flg) { 50116d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxWt(jac->hsolver,tmpdbl);CHKERRQ(ierr); 50216d9e3a6SLisandro Dalcin jac->relaxweight = tmpdbl; 50316d9e3a6SLisandro Dalcin } 50416d9e3a6SLisandro Dalcin 50516d9e3a6SLisandro Dalcin n=2; 50616d9e3a6SLisandro Dalcin twodbl[0] = twodbl[1] = 1.0; 50716d9e3a6SLisandro 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); 50816d9e3a6SLisandro Dalcin if (flg) { 50916d9e3a6SLisandro Dalcin if (n == 2) { 51016d9e3a6SLisandro Dalcin indx = (int)PetscAbsReal(twodbl[1]); 51116d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetLevelRelaxWt(jac->hsolver,twodbl[0],indx);CHKERRQ(ierr); 51216d9e3a6SLisandro Dalcin } else { 51316d9e3a6SLisandro 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); 51416d9e3a6SLisandro Dalcin } 51516d9e3a6SLisandro Dalcin } 51616d9e3a6SLisandro Dalcin 51716d9e3a6SLisandro Dalcin /* Outer relaxation Weight */ 51816d9e3a6SLisandro 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); 51916d9e3a6SLisandro Dalcin if (flg) { 52016d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetOuterWt( jac->hsolver, tmpdbl);CHKERRQ(ierr); 52116d9e3a6SLisandro Dalcin jac->outerrelaxweight = tmpdbl; 52216d9e3a6SLisandro Dalcin } 52316d9e3a6SLisandro Dalcin 52416d9e3a6SLisandro Dalcin n=2; 52516d9e3a6SLisandro Dalcin twodbl[0] = twodbl[1] = 1.0; 52616d9e3a6SLisandro 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); 52716d9e3a6SLisandro Dalcin if (flg) { 52816d9e3a6SLisandro Dalcin if (n == 2) { 52916d9e3a6SLisandro Dalcin indx = (int)PetscAbsReal(twodbl[1]); 53016d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetLevelOuterWt( jac->hsolver, twodbl[0], indx);CHKERRQ(ierr); 53116d9e3a6SLisandro Dalcin } else { 53216d9e3a6SLisandro 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); 53316d9e3a6SLisandro Dalcin } 53416d9e3a6SLisandro Dalcin } 53516d9e3a6SLisandro Dalcin 53616d9e3a6SLisandro Dalcin /* the Relax Order */ 53716d9e3a6SLisandro Dalcin /* ierr = PetscOptionsName("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", &flg);CHKERRQ(ierr); */ 53816d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr); 53916d9e3a6SLisandro Dalcin 54016d9e3a6SLisandro Dalcin if (flg) { 54116d9e3a6SLisandro Dalcin jac->relaxorder = 0; 54216d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr); 54316d9e3a6SLisandro Dalcin } 54416d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,2,HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr); 54516d9e3a6SLisandro Dalcin if (flg) { 54616d9e3a6SLisandro Dalcin jac->measuretype = indx; 54716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr); 54816d9e3a6SLisandro Dalcin } 5490f1074feSSatish Balay /* update list length 3/07 */ 5500f1074feSSatish Balay ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,11,HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr); 55116d9e3a6SLisandro Dalcin if (flg) { 55216d9e3a6SLisandro Dalcin jac->coarsentype = indx; 55316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr); 55416d9e3a6SLisandro Dalcin } 5550f1074feSSatish Balay 5560f1074feSSatish Balay /* new 3/07 */ 5570f1074feSSatish Balay ierr = PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,14,HYPREBoomerAMGInterpType[0],&indx,&flg);CHKERRQ(ierr); 5580f1074feSSatish Balay if (flg) { 5590f1074feSSatish Balay jac->interptype = indx; 5600f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);CHKERRQ(ierr); 5610f1074feSSatish Balay } 5620f1074feSSatish Balay 5630f1074feSSatish Balay 56416d9e3a6SLisandro Dalcin ierr = PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);CHKERRQ(ierr); 56516d9e3a6SLisandro Dalcin if (flg) { 56616d9e3a6SLisandro Dalcin int level=3; 56716d9e3a6SLisandro Dalcin jac->printstatistics = PETSC_TRUE; 56816d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,PETSC_NULL);CHKERRQ(ierr); 56916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetPrintLevel(jac->hsolver,level);CHKERRQ(ierr); 570*2ae77aedSBarry Smith } 571*2ae77aedSBarry Smith 572*2ae77aedSBarry Smith ierr = PetscOptionsName("-pc_hypre_boomeramg_print_debug","Print debug information","None",&flg);CHKERRQ(ierr); 573*2ae77aedSBarry Smith if (flg) { 574*2ae77aedSBarry Smith int level=3; 575*2ae77aedSBarry Smith jac->printstatistics = PETSC_TRUE; 576*2ae77aedSBarry Smith ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,PETSC_NULL);CHKERRQ(ierr); 57716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetDebugFlag(jac->hsolver,level);CHKERRQ(ierr); 57816d9e3a6SLisandro Dalcin } 5798f87f92bSBarry Smith 5808f87f92bSBarry Smith ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr); 5818f87f92bSBarry Smith if (flg && tmp_truth) { 5828f87f92bSBarry Smith jac->nodal_coarsen = 1; 5838f87f92bSBarry Smith ierr = HYPRE_BoomerAMGSetNodal(jac->hsolver,1);CHKERRQ(ierr); 5848f87f92bSBarry Smith } 5858f87f92bSBarry Smith 5868f87f92bSBarry Smith ierr = PetscOptionsTruth( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr); 5878f87f92bSBarry Smith if (flg && tmp_truth) { 5888f87f92bSBarry Smith PetscInt tmp_int; 5898f87f92bSBarry Smith ierr = PetscOptionsInt( "-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);CHKERRQ(ierr); 5908f87f92bSBarry Smith if (flg) jac->nodal_relax_levels = tmp_int; 5918f87f92bSBarry Smith ierr = HYPRE_BoomerAMGSetSmoothType(jac->hsolver,6);CHKERRQ(ierr); 5928f87f92bSBarry Smith ierr = HYPRE_BoomerAMGSetDomainType(jac->hsolver,1);CHKERRQ(ierr); 5938f87f92bSBarry Smith ierr = HYPRE_BoomerAMGSetOverlap(jac->hsolver,0);CHKERRQ(ierr); 5948f87f92bSBarry Smith ierr = HYPRE_BoomerAMGSetSmoothNumLevels(jac->hsolver,jac->nodal_relax_levels);CHKERRQ(ierr); 5958f87f92bSBarry Smith } 5968f87f92bSBarry Smith 59716d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 59816d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 59916d9e3a6SLisandro Dalcin } 60016d9e3a6SLisandro Dalcin 60116d9e3a6SLisandro Dalcin #undef __FUNCT__ 60216d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG" 60316d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its) 60416d9e3a6SLisandro Dalcin { 60516d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 60616d9e3a6SLisandro Dalcin PetscErrorCode ierr; 60716d9e3a6SLisandro Dalcin 60816d9e3a6SLisandro Dalcin PetscFunctionBegin; 60916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,its*jac->maxiter);CHKERRQ(ierr); 61016d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,PetscMin(rtol,jac->tol));CHKERRQ(ierr); 61116d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_TRUE; 61216d9e3a6SLisandro Dalcin ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr); 61316d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_FALSE; 61416d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 61516d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 61616d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 61716d9e3a6SLisandro Dalcin } 61816d9e3a6SLisandro Dalcin 61916d9e3a6SLisandro Dalcin 62016d9e3a6SLisandro Dalcin #undef __FUNCT__ 62116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_BoomerAMG" 62216d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer) 62316d9e3a6SLisandro Dalcin { 62416d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 62516d9e3a6SLisandro Dalcin PetscErrorCode ierr; 62616d9e3a6SLisandro Dalcin PetscTruth iascii; 62716d9e3a6SLisandro Dalcin 62816d9e3a6SLisandro Dalcin PetscFunctionBegin; 62916d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 63016d9e3a6SLisandro Dalcin if (iascii) { 63116d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr); 63216d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr); 63316d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr); 63416d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr); 63516d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Convergence tolerance PER hypre call %G\n",jac->tol);CHKERRQ(ierr); 63616d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Threshold for strong coupling %G\n",jac->strongthreshold);CHKERRQ(ierr); 6370f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation truncation factor %G\n",jac->truncfactor);CHKERRQ(ierr); 6380f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);CHKERRQ(ierr); 6390f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);CHKERRQ(ierr); 6400f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);CHKERRQ(ierr); 6410f1074feSSatish Balay 64216d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum row sums %G\n",jac->maxrowsum);CHKERRQ(ierr); 64316d9e3a6SLisandro Dalcin 6440f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps down %d\n",jac->gridsweeps[0]);CHKERRQ(ierr); 6450f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps up %d\n",jac->gridsweeps[1]);CHKERRQ(ierr); 6460f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps on coarse %d\n",jac->gridsweeps[2]);CHKERRQ(ierr); 64716d9e3a6SLisandro Dalcin 6480f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax down %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr); 6490f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax up %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr); 6500f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax on coarse %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr); 65116d9e3a6SLisandro Dalcin 65216d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax weight (all) %G\n",jac->relaxweight);CHKERRQ(ierr); 65316d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Outer relax weight (all) %G\n",jac->outerrelaxweight);CHKERRQ(ierr); 65416d9e3a6SLisandro Dalcin 65516d9e3a6SLisandro Dalcin if (jac->relaxorder) { 65616d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr); 65716d9e3a6SLisandro Dalcin } else { 65816d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr); 65916d9e3a6SLisandro Dalcin } 66016d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Measure type %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr); 66116d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Coarsen type %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr); 6620f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation type %s\n",HYPREBoomerAMGInterpType[jac->interptype]);CHKERRQ(ierr); 6638f87f92bSBarry Smith if (jac->nodal_coarsen) { 6648f87f92bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal())\n");CHKERRQ(ierr); 6658f87f92bSBarry Smith } 6668f87f92bSBarry Smith if (jac->nodal_relax) { 6678f87f92bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);CHKERRQ(ierr); 6688f87f92bSBarry Smith } 66916d9e3a6SLisandro Dalcin } 67016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 67116d9e3a6SLisandro Dalcin } 67216d9e3a6SLisandro Dalcin 67316d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 67416d9e3a6SLisandro Dalcin #undef __FUNCT__ 67516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails" 67616d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PC pc) 67716d9e3a6SLisandro Dalcin { 67816d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 67916d9e3a6SLisandro Dalcin PetscErrorCode ierr; 68016d9e3a6SLisandro Dalcin int indx; 68116d9e3a6SLisandro Dalcin PetscTruth flag; 68216d9e3a6SLisandro Dalcin const char *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"}; 68316d9e3a6SLisandro Dalcin 68416d9e3a6SLisandro Dalcin PetscFunctionBegin; 68516d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE ParaSails Options");CHKERRQ(ierr); 68616d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr); 68716d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr); 68816d9e3a6SLisandro Dalcin if (flag) { 68916d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr); 69016d9e3a6SLisandro Dalcin } 69116d9e3a6SLisandro Dalcin 69216d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr); 69316d9e3a6SLisandro Dalcin if (flag) { 69416d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr); 69516d9e3a6SLisandro Dalcin } 69616d9e3a6SLisandro Dalcin 69716d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr); 69816d9e3a6SLisandro Dalcin if (flag) { 69916d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr); 70016d9e3a6SLisandro Dalcin } 70116d9e3a6SLisandro Dalcin 70216d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_parasails_logging","Print info to screen","None",(PetscTruth)jac->logging,(PetscTruth*)&jac->logging,&flag);CHKERRQ(ierr); 70316d9e3a6SLisandro Dalcin if (flag) { 70416d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr); 70516d9e3a6SLisandro Dalcin } 70616d9e3a6SLisandro Dalcin 70716d9e3a6SLisandro Dalcin ierr = PetscOptionsTruth("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscTruth)jac->ruse,(PetscTruth*)&jac->ruse,&flag);CHKERRQ(ierr); 70816d9e3a6SLisandro Dalcin if (flag) { 70916d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr); 71016d9e3a6SLisandro Dalcin } 71116d9e3a6SLisandro Dalcin 71216d9e3a6SLisandro Dalcin ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,3,symtlist[0],&indx,&flag);CHKERRQ(ierr); 71316d9e3a6SLisandro Dalcin if (flag) { 71416d9e3a6SLisandro Dalcin jac->symt = indx; 71516d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr); 71616d9e3a6SLisandro Dalcin } 71716d9e3a6SLisandro Dalcin 71816d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 71916d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 72016d9e3a6SLisandro Dalcin } 72116d9e3a6SLisandro Dalcin 72216d9e3a6SLisandro Dalcin #undef __FUNCT__ 72316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_ParaSails" 72416d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer) 72516d9e3a6SLisandro Dalcin { 72616d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 72716d9e3a6SLisandro Dalcin PetscErrorCode ierr; 72816d9e3a6SLisandro Dalcin PetscTruth iascii; 72916d9e3a6SLisandro Dalcin const char *symt = 0;; 73016d9e3a6SLisandro Dalcin 73116d9e3a6SLisandro Dalcin PetscFunctionBegin; 73216d9e3a6SLisandro Dalcin ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 73316d9e3a6SLisandro Dalcin if (iascii) { 73416d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails preconditioning\n");CHKERRQ(ierr); 73516d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr); 73616d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: threshold %G\n",jac->threshhold);CHKERRQ(ierr); 73716d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: filter %G\n",jac->filter);CHKERRQ(ierr); 73816d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: load balance %G\n",jac->loadbal);CHKERRQ(ierr); 73916d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: reuse nonzero structure %s\n",PetscTruths[jac->ruse]);CHKERRQ(ierr); 74016d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: print info to screen %s\n",PetscTruths[jac->logging]);CHKERRQ(ierr); 74116d9e3a6SLisandro Dalcin if (!jac->symt) { 74216d9e3a6SLisandro Dalcin symt = "nonsymmetric matrix and preconditioner"; 74316d9e3a6SLisandro Dalcin } else if (jac->symt == 1) { 74416d9e3a6SLisandro Dalcin symt = "SPD matrix and preconditioner"; 74516d9e3a6SLisandro Dalcin } else if (jac->symt == 2) { 74616d9e3a6SLisandro Dalcin symt = "nonsymmetric matrix but SPD preconditioner"; 74716d9e3a6SLisandro Dalcin } else { 74816d9e3a6SLisandro Dalcin SETERRQ1(PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt); 74916d9e3a6SLisandro Dalcin } 75016d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr); 75116d9e3a6SLisandro Dalcin } 75216d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 75316d9e3a6SLisandro Dalcin } 75416d9e3a6SLisandro Dalcin /* ---------------------------------------------------------------------------------*/ 75516d9e3a6SLisandro Dalcin 75616d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 75716d9e3a6SLisandro Dalcin #undef __FUNCT__ 75816d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType_HYPRE" 75916d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType_HYPRE(PC pc,const char *name[]) 76016d9e3a6SLisandro Dalcin { 76116d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 76216d9e3a6SLisandro Dalcin 76316d9e3a6SLisandro Dalcin PetscFunctionBegin; 76416d9e3a6SLisandro Dalcin *name = jac->hypre_type; 76516d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 76616d9e3a6SLisandro Dalcin } 76716d9e3a6SLisandro Dalcin EXTERN_C_END 76816d9e3a6SLisandro Dalcin 76916d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 77016d9e3a6SLisandro Dalcin #undef __FUNCT__ 77116d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType_HYPRE" 77216d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType_HYPRE(PC pc,const char name[]) 77316d9e3a6SLisandro Dalcin { 77416d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 77516d9e3a6SLisandro Dalcin PetscErrorCode ierr; 77616d9e3a6SLisandro Dalcin PetscTruth flag; 77716d9e3a6SLisandro Dalcin 77816d9e3a6SLisandro Dalcin PetscFunctionBegin; 77916d9e3a6SLisandro Dalcin if (jac->hypre_type) { 78016d9e3a6SLisandro Dalcin ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr); 78116d9e3a6SLisandro Dalcin if (!flag) SETERRQ(PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set"); 78216d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 78316d9e3a6SLisandro Dalcin } else { 78416d9e3a6SLisandro Dalcin ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr); 78516d9e3a6SLisandro Dalcin } 78616d9e3a6SLisandro Dalcin 78716d9e3a6SLisandro Dalcin jac->maxiter = PETSC_DEFAULT; 78816d9e3a6SLisandro Dalcin jac->tol = PETSC_DEFAULT; 78916d9e3a6SLisandro Dalcin jac->printstatistics = PetscLogPrintInfo; 79016d9e3a6SLisandro Dalcin 79116d9e3a6SLisandro Dalcin ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr); 79216d9e3a6SLisandro Dalcin if (flag) { 79316d9e3a6SLisandro Dalcin ierr = HYPRE_ParCSRPilutCreate(jac->comm_hypre,&jac->hsolver); 79416d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut; 79516d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_Pilut; 79616d9e3a6SLisandro Dalcin jac->destroy = HYPRE_ParCSRPilutDestroy; 79716d9e3a6SLisandro Dalcin jac->setup = HYPRE_ParCSRPilutSetup; 79816d9e3a6SLisandro Dalcin jac->solve = HYPRE_ParCSRPilutSolve; 79916d9e3a6SLisandro Dalcin jac->factorrowsize = PETSC_DEFAULT; 80016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 80116d9e3a6SLisandro Dalcin } 80216d9e3a6SLisandro Dalcin ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr); 80316d9e3a6SLisandro Dalcin if (flag) { 80416d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsCreate(jac->comm_hypre,&jac->hsolver); 80516d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails; 80616d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_ParaSails; 80716d9e3a6SLisandro Dalcin jac->destroy = HYPRE_ParaSailsDestroy; 80816d9e3a6SLisandro Dalcin jac->setup = HYPRE_ParaSailsSetup; 80916d9e3a6SLisandro Dalcin jac->solve = HYPRE_ParaSailsSolve; 81016d9e3a6SLisandro Dalcin /* initialize */ 81116d9e3a6SLisandro Dalcin jac->nlevels = 1; 81216d9e3a6SLisandro Dalcin jac->threshhold = .1; 81316d9e3a6SLisandro Dalcin jac->filter = .1; 81416d9e3a6SLisandro Dalcin jac->loadbal = 0; 81516d9e3a6SLisandro Dalcin if (PetscLogPrintInfo) { 81616d9e3a6SLisandro Dalcin jac->logging = (int) PETSC_TRUE; 81716d9e3a6SLisandro Dalcin } else { 81816d9e3a6SLisandro Dalcin jac->logging = (int) PETSC_FALSE; 81916d9e3a6SLisandro Dalcin } 82016d9e3a6SLisandro Dalcin jac->ruse = (int) PETSC_FALSE; 82116d9e3a6SLisandro Dalcin jac->symt = 0; 82216d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetParams(jac->hsolver,jac->threshhold,jac->nlevels);CHKERRQ(ierr); 82316d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetFilter(jac->hsolver,jac->filter);CHKERRQ(ierr); 82416d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLoadbal(jac->hsolver,jac->loadbal);CHKERRQ(ierr); 82516d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetLogging(jac->hsolver,jac->logging);CHKERRQ(ierr); 82616d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetReuse(jac->hsolver,jac->ruse);CHKERRQ(ierr); 82716d9e3a6SLisandro Dalcin ierr = HYPRE_ParaSailsSetSym(jac->hsolver,jac->symt);CHKERRQ(ierr); 82816d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 82916d9e3a6SLisandro Dalcin } 83016d9e3a6SLisandro Dalcin ierr = PetscStrcmp("euclid",jac->hypre_type,&flag);CHKERRQ(ierr); 83116d9e3a6SLisandro Dalcin if (flag) { 83216d9e3a6SLisandro Dalcin ierr = HYPRE_EuclidCreate(jac->comm_hypre,&jac->hsolver); 83316d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Euclid; 83416d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_Euclid; 83516d9e3a6SLisandro Dalcin jac->destroy = HYPRE_EuclidDestroy; 83616d9e3a6SLisandro Dalcin jac->setup = HYPRE_EuclidSetup; 83716d9e3a6SLisandro Dalcin jac->solve = HYPRE_EuclidSolve; 83816d9e3a6SLisandro Dalcin /* initialization */ 83916d9e3a6SLisandro Dalcin jac->bjilu = PETSC_FALSE; 84016d9e3a6SLisandro Dalcin jac->levels = 1; 84116d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 84216d9e3a6SLisandro Dalcin } 84316d9e3a6SLisandro Dalcin ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr); 84416d9e3a6SLisandro Dalcin if (flag) { 84516d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGCreate(&jac->hsolver); 84616d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_BoomerAMG; 84716d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_BoomerAMG; 84816d9e3a6SLisandro Dalcin pc->ops->applytranspose = PCApplyTranspose_HYPRE_BoomerAMG; 84916d9e3a6SLisandro Dalcin pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG; 85016d9e3a6SLisandro Dalcin jac->destroy = HYPRE_BoomerAMGDestroy; 85116d9e3a6SLisandro Dalcin jac->setup = HYPRE_BoomerAMGSetup; 85216d9e3a6SLisandro Dalcin jac->solve = HYPRE_BoomerAMGSolve; 85316d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_FALSE; 85416d9e3a6SLisandro Dalcin /* these defaults match the hypre defaults */ 85516d9e3a6SLisandro Dalcin jac->cycletype = 1; 85616d9e3a6SLisandro Dalcin jac->maxlevels = 25; 85716d9e3a6SLisandro Dalcin jac->maxiter = 1; 8588f87f92bSBarry Smith jac->tol = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */ 85916d9e3a6SLisandro Dalcin jac->truncfactor = 0.0; 86016d9e3a6SLisandro Dalcin jac->strongthreshold = .25; 86116d9e3a6SLisandro Dalcin jac->maxrowsum = .9; 86216d9e3a6SLisandro Dalcin jac->coarsentype = 6; 86316d9e3a6SLisandro Dalcin jac->measuretype = 0; 8640f1074feSSatish Balay jac->gridsweeps[0] = jac->gridsweeps[1] = jac->gridsweeps[2] = 1; 8658f87f92bSBarry Smith jac->relaxtype[0] = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */ 8660f1074feSSatish Balay jac->relaxtype[2] = 9; /*G.E. */ 86716d9e3a6SLisandro Dalcin jac->relaxweight = 1.0; 86816d9e3a6SLisandro Dalcin jac->outerrelaxweight = 1.0; 86916d9e3a6SLisandro Dalcin jac->relaxorder = 1; 8700f1074feSSatish Balay jac->interptype = 0; 8710f1074feSSatish Balay jac->agg_nl = 0; 8720f1074feSSatish Balay jac->pmax = 0; 8730f1074feSSatish Balay jac->truncfactor = 0.0; 8740f1074feSSatish Balay jac->agg_num_paths = 1; 8758f87f92bSBarry Smith 8768f87f92bSBarry Smith jac->nodal_coarsen = 0; 8778f87f92bSBarry Smith jac->nodal_relax = PETSC_FALSE; 8788f87f92bSBarry Smith jac->nodal_relax_levels = 1; 87916d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCycleType(jac->hsolver,jac->cycletype);CHKERRQ(ierr); 88016d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxLevels(jac->hsolver,jac->maxlevels);CHKERRQ(ierr); 88116d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxIter(jac->hsolver,jac->maxiter);CHKERRQ(ierr); 88216d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTol(jac->hsolver,jac->tol);CHKERRQ(ierr); 88316d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetTruncFactor(jac->hsolver,jac->truncfactor);CHKERRQ(ierr); 88416d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetStrongThreshold(jac->hsolver,jac->strongthreshold);CHKERRQ(ierr); 88516d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMaxRowSum(jac->hsolver,jac->maxrowsum);CHKERRQ(ierr); 88616d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetCoarsenType(jac->hsolver,jac->coarsentype);CHKERRQ(ierr); 88716d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetMeasureType(jac->hsolver,jac->measuretype);CHKERRQ(ierr); 88816d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGSetRelaxOrder(jac->hsolver, jac->relaxorder);CHKERRQ(ierr); 8890f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetInterpType(jac->hsolver,jac->interptype);CHKERRQ(ierr); 8900f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetAggNumLevels(jac->hsolver,jac->agg_nl); 8910f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetPMaxElmts(jac->hsolver,jac->pmax);CHKERRQ(ierr); 8920f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetNumPaths(jac->hsolver,jac->agg_num_paths);CHKERRQ(ierr); 8930f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetRelaxType(jac->hsolver, jac->relaxtype[0]); /*defaults coarse to 9*/ 8940f1074feSSatish Balay ierr = HYPRE_BoomerAMGSetNumSweeps(jac->hsolver, jac->gridsweeps[0]); /*defaults coarse to 1 */ 89516d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 89616d9e3a6SLisandro Dalcin } 89716d9e3a6SLisandro Dalcin ierr = PetscStrfree(jac->hypre_type);CHKERRQ(ierr); 89816d9e3a6SLisandro Dalcin jac->hypre_type = PETSC_NULL; 89916d9e3a6SLisandro Dalcin SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, euclid, boomeramg",name); 90016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 90116d9e3a6SLisandro Dalcin } 90216d9e3a6SLisandro Dalcin EXTERN_C_END 90316d9e3a6SLisandro Dalcin 90416d9e3a6SLisandro Dalcin /* 90516d9e3a6SLisandro Dalcin It only gets here if the HYPRE type has not been set before the call to 90616d9e3a6SLisandro Dalcin ...SetFromOptions() which actually is most of the time 90716d9e3a6SLisandro Dalcin */ 90816d9e3a6SLisandro Dalcin #undef __FUNCT__ 90916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE" 91016d9e3a6SLisandro Dalcin static PetscErrorCode PCSetFromOptions_HYPRE(PC pc) 91116d9e3a6SLisandro Dalcin { 91216d9e3a6SLisandro Dalcin PetscErrorCode ierr; 91316d9e3a6SLisandro Dalcin int indx; 91416d9e3a6SLisandro Dalcin const char *type[] = {"pilut","parasails","boomeramg","euclid"}; 91516d9e3a6SLisandro Dalcin PetscTruth flg; 91616d9e3a6SLisandro Dalcin 91716d9e3a6SLisandro Dalcin PetscFunctionBegin; 91816d9e3a6SLisandro Dalcin ierr = PetscOptionsHead("HYPRE preconditioner options");CHKERRQ(ierr); 91902a17cd4SBarry Smith ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,4,"boomeramg",&indx,&flg);CHKERRQ(ierr); 92016d9e3a6SLisandro Dalcin if (flg) { 92116d9e3a6SLisandro Dalcin ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr); 92202a17cd4SBarry Smith } else { 92302a17cd4SBarry Smith ierr = PCHYPRESetType_HYPRE(pc,"boomeramg");CHKERRQ(ierr); 92416d9e3a6SLisandro Dalcin } 92516d9e3a6SLisandro Dalcin if (pc->ops->setfromoptions) { 92616d9e3a6SLisandro Dalcin ierr = pc->ops->setfromoptions(pc);CHKERRQ(ierr); 92716d9e3a6SLisandro Dalcin } 92816d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 92916d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 93016d9e3a6SLisandro Dalcin } 93116d9e3a6SLisandro Dalcin 93216d9e3a6SLisandro Dalcin #undef __FUNCT__ 93316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType" 93416d9e3a6SLisandro Dalcin /*@C 93516d9e3a6SLisandro Dalcin PCHYPRESetType - Sets which hypre preconditioner you wish to use 93616d9e3a6SLisandro Dalcin 93716d9e3a6SLisandro Dalcin Input Parameters: 93816d9e3a6SLisandro Dalcin + pc - the preconditioner context 93916d9e3a6SLisandro Dalcin - name - either pilut, parasails, boomeramg, euclid 94016d9e3a6SLisandro Dalcin 94116d9e3a6SLisandro Dalcin Options Database Keys: 94216d9e3a6SLisandro Dalcin -pc_hypre_type - One of pilut, parasails, boomeramg, euclid 94316d9e3a6SLisandro Dalcin 94416d9e3a6SLisandro Dalcin Level: intermediate 94516d9e3a6SLisandro Dalcin 94616d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 94716d9e3a6SLisandro Dalcin PCHYPRE 94816d9e3a6SLisandro Dalcin 94916d9e3a6SLisandro Dalcin @*/ 95016d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPRESetType(PC pc,const char name[]) 95116d9e3a6SLisandro Dalcin { 95216d9e3a6SLisandro Dalcin PetscErrorCode ierr,(*f)(PC,const char[]); 95316d9e3a6SLisandro Dalcin 95416d9e3a6SLisandro Dalcin PetscFunctionBegin; 95516d9e3a6SLisandro Dalcin PetscValidHeaderSpecific(pc,PC_COOKIE,1); 95616d9e3a6SLisandro Dalcin PetscValidCharPointer(name,2); 95716d9e3a6SLisandro Dalcin ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPRESetType_C",(void (**)(void))&f);CHKERRQ(ierr); 95816d9e3a6SLisandro Dalcin if (f) { 95916d9e3a6SLisandro Dalcin ierr = (*f)(pc,name);CHKERRQ(ierr); 96016d9e3a6SLisandro Dalcin } 96116d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 96216d9e3a6SLisandro Dalcin } 96316d9e3a6SLisandro Dalcin 96416d9e3a6SLisandro Dalcin #undef __FUNCT__ 96516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType" 96616d9e3a6SLisandro Dalcin /*@C 96716d9e3a6SLisandro Dalcin PCHYPREGetType - Gets which hypre preconditioner you are using 96816d9e3a6SLisandro Dalcin 96916d9e3a6SLisandro Dalcin Input Parameter: 97016d9e3a6SLisandro Dalcin . pc - the preconditioner context 97116d9e3a6SLisandro Dalcin 97216d9e3a6SLisandro Dalcin Output Parameter: 97316d9e3a6SLisandro Dalcin . name - either pilut, parasails, boomeramg, euclid 97416d9e3a6SLisandro Dalcin 97516d9e3a6SLisandro Dalcin Level: intermediate 97616d9e3a6SLisandro Dalcin 97716d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC, 97816d9e3a6SLisandro Dalcin PCHYPRE 97916d9e3a6SLisandro Dalcin 98016d9e3a6SLisandro Dalcin @*/ 98116d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCHYPREGetType(PC pc,const char *name[]) 98216d9e3a6SLisandro Dalcin { 98316d9e3a6SLisandro Dalcin PetscErrorCode ierr,(*f)(PC,const char*[]); 98416d9e3a6SLisandro Dalcin 98516d9e3a6SLisandro Dalcin PetscFunctionBegin; 98616d9e3a6SLisandro Dalcin PetscValidHeaderSpecific(pc,PC_COOKIE,1); 98716d9e3a6SLisandro Dalcin PetscValidPointer(name,2); 98816d9e3a6SLisandro Dalcin ierr = PetscObjectQueryFunction((PetscObject)pc,"PCHYPREGetType_C",(void (**)(void))&f);CHKERRQ(ierr); 98916d9e3a6SLisandro Dalcin if (f) { 99016d9e3a6SLisandro Dalcin ierr = (*f)(pc,name);CHKERRQ(ierr); 99116d9e3a6SLisandro Dalcin } 99216d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 99316d9e3a6SLisandro Dalcin } 99416d9e3a6SLisandro Dalcin 99516d9e3a6SLisandro Dalcin /*MC 99616d9e3a6SLisandro Dalcin PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre 99716d9e3a6SLisandro Dalcin 99816d9e3a6SLisandro Dalcin Options Database Keys: 99916d9e3a6SLisandro Dalcin + -pc_hypre_type - One of pilut, parasails, boomeramg, euclid 100016d9e3a6SLisandro Dalcin - Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX 100116d9e3a6SLisandro Dalcin preconditioner 100216d9e3a6SLisandro Dalcin 100316d9e3a6SLisandro Dalcin Level: intermediate 100416d9e3a6SLisandro Dalcin 100516d9e3a6SLisandro Dalcin Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()), 100616d9e3a6SLisandro Dalcin the many hypre options can ONLY be set via the options database (e.g. the command line 100716d9e3a6SLisandro Dalcin or with PetscOptionsSetValue(), there are no functions to set them) 100816d9e3a6SLisandro Dalcin 100916d9e3a6SLisandro Dalcin The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations 10100f1074feSSatish Balay (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if 10110f1074feSSatish Balay -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner 10120f1074feSSatish Balay (-pc_hypre_boomeramg_rtol should be set to 0.0 - the default - to strictly use a fixed number of 10138f87f92bSBarry Smith iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations 10140f1074feSSatish Balay and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10 10150f1074feSSatish Balay then AT MOST twenty V-cycles of boomeramg will be called. 101616d9e3a6SLisandro Dalcin 10170f1074feSSatish Balay Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation 10180f1074feSSatish Balay (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry. 10190f1074feSSatish Balay Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi. 102016d9e3a6SLisandro Dalcin If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly 102116d9e3a6SLisandro Dalcin and use -ksp_max_it to control the number of V-cycles. 102216d9e3a6SLisandro Dalcin (see the PETSc FAQ.html at the PETSc website under the Documentation tab). 102316d9e3a6SLisandro Dalcin 102416d9e3a6SLisandro Dalcin 2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option 102516d9e3a6SLisandro Dalcin -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L. 102616d9e3a6SLisandro Dalcin 102716d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 102816d9e3a6SLisandro Dalcin PCHYPRESetType() 102916d9e3a6SLisandro Dalcin 103016d9e3a6SLisandro Dalcin M*/ 103116d9e3a6SLisandro Dalcin 103216d9e3a6SLisandro Dalcin EXTERN_C_BEGIN 103316d9e3a6SLisandro Dalcin #undef __FUNCT__ 103416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCCreate_HYPRE" 103516d9e3a6SLisandro Dalcin PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_HYPRE(PC pc) 103616d9e3a6SLisandro Dalcin { 103716d9e3a6SLisandro Dalcin PC_HYPRE *jac; 103816d9e3a6SLisandro Dalcin PetscErrorCode ierr; 103916d9e3a6SLisandro Dalcin 104016d9e3a6SLisandro Dalcin PetscFunctionBegin; 104138f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_HYPRE,&jac);CHKERRQ(ierr); 104216d9e3a6SLisandro Dalcin pc->data = jac; 104316d9e3a6SLisandro Dalcin pc->ops->destroy = PCDestroy_HYPRE; 104416d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE; 104516d9e3a6SLisandro Dalcin pc->ops->setup = PCSetUp_HYPRE; 104616d9e3a6SLisandro Dalcin pc->ops->apply = PCApply_HYPRE; 104716d9e3a6SLisandro Dalcin jac->comm_hypre = MPI_COMM_NULL; 104816d9e3a6SLisandro Dalcin jac->hypre_type = PETSC_NULL; 104916d9e3a6SLisandro Dalcin /* duplicate communicator for hypre */ 10507adad957SLisandro Dalcin ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(jac->comm_hypre));CHKERRQ(ierr); 105116d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPRESetType_C","PCHYPRESetType_HYPRE",PCHYPRESetType_HYPRE);CHKERRQ(ierr); 105216d9e3a6SLisandro Dalcin ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCHYPREGetType_C","PCHYPREGetType_HYPRE",PCHYPREGetType_HYPRE);CHKERRQ(ierr); 105316d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 105416d9e3a6SLisandro Dalcin } 105516d9e3a6SLisandro Dalcin EXTERN_C_END 1056