116d9e3a6SLisandro Dalcin 216d9e3a6SLisandro Dalcin /* 316d9e3a6SLisandro Dalcin Provides an interface to the LLNL package hypre 416d9e3a6SLisandro Dalcin */ 50f1074feSSatish Balay 60f1074feSSatish Balay /* Must use hypre 2.0.0 or more recent. */ 70f1074feSSatish Balay 8b45d2f2cSJed Brown #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/ 9c6db04a5SJed Brown #include <../src/dm/impls/da/hypre/mhyp.h> 1016d9e3a6SLisandro Dalcin 11dff31646SBarry Smith static PetscBool cite = PETSC_FALSE; 121f817a21SBarry Smith static const char hypreCitation[] = "@manual{hypre-web-page,\n title = {{\\sl hypre}: High Performance Preconditioners},\n organization = {Lawrence Livermore National Laboratory},\n note = {\\url{http://www.llnl.gov/CASC/hypre/}}\n}\n"; 131f817a21SBarry Smith 1416d9e3a6SLisandro Dalcin /* 1516d9e3a6SLisandro Dalcin Private context (data structure) for the preconditioner. 1616d9e3a6SLisandro Dalcin */ 1716d9e3a6SLisandro Dalcin typedef struct { 1816d9e3a6SLisandro Dalcin HYPRE_Solver hsolver; 1916d9e3a6SLisandro Dalcin HYPRE_IJMatrix ij; 2016d9e3a6SLisandro Dalcin HYPRE_IJVector b,x; 2116d9e3a6SLisandro Dalcin 224ddd07fcSJed Brown HYPRE_Int (*destroy)(HYPRE_Solver); 234ddd07fcSJed Brown HYPRE_Int (*solve)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector); 244ddd07fcSJed Brown HYPRE_Int (*setup)(HYPRE_Solver,HYPRE_ParCSRMatrix,HYPRE_ParVector,HYPRE_ParVector); 2516d9e3a6SLisandro Dalcin 2616d9e3a6SLisandro Dalcin MPI_Comm comm_hypre; 2716d9e3a6SLisandro Dalcin char *hypre_type; 2816d9e3a6SLisandro Dalcin 2916d9e3a6SLisandro Dalcin /* options for Pilut and BoomerAMG*/ 304ddd07fcSJed Brown PetscInt maxiter; 3116d9e3a6SLisandro Dalcin double tol; 3216d9e3a6SLisandro Dalcin 3316d9e3a6SLisandro Dalcin /* options for Pilut */ 344ddd07fcSJed Brown PetscInt factorrowsize; 3516d9e3a6SLisandro Dalcin 3616d9e3a6SLisandro Dalcin /* options for ParaSails */ 374ddd07fcSJed Brown PetscInt nlevels; 3816d9e3a6SLisandro Dalcin double threshhold; 3916d9e3a6SLisandro Dalcin double filter; 404ddd07fcSJed Brown PetscInt sym; 4116d9e3a6SLisandro Dalcin double loadbal; 424ddd07fcSJed Brown PetscInt logging; 434ddd07fcSJed Brown PetscInt ruse; 444ddd07fcSJed Brown PetscInt symt; 4516d9e3a6SLisandro Dalcin 46*22b6d1caSBarry Smith /* options for BoomerAMG */ 47ace3abfcSBarry Smith PetscBool printstatistics; 4816d9e3a6SLisandro Dalcin 4916d9e3a6SLisandro Dalcin /* options for BoomerAMG */ 504ddd07fcSJed Brown PetscInt cycletype; 514ddd07fcSJed Brown PetscInt maxlevels; 5216d9e3a6SLisandro Dalcin double strongthreshold; 5316d9e3a6SLisandro Dalcin double maxrowsum; 544ddd07fcSJed Brown PetscInt gridsweeps[3]; 554ddd07fcSJed Brown PetscInt coarsentype; 564ddd07fcSJed Brown PetscInt measuretype; 574ddd07fcSJed Brown PetscInt relaxtype[3]; 5816d9e3a6SLisandro Dalcin double relaxweight; 5916d9e3a6SLisandro Dalcin double outerrelaxweight; 604ddd07fcSJed Brown PetscInt relaxorder; 6116d9e3a6SLisandro Dalcin double truncfactor; 62ace3abfcSBarry Smith PetscBool applyrichardson; 634ddd07fcSJed Brown PetscInt pmax; 644ddd07fcSJed Brown PetscInt interptype; 654ddd07fcSJed Brown PetscInt agg_nl; 664ddd07fcSJed Brown PetscInt agg_num_paths; 674ddd07fcSJed Brown PetscInt nodal_coarsen; 68ace3abfcSBarry Smith PetscBool nodal_relax; 694ddd07fcSJed Brown PetscInt nodal_relax_levels; 7016d9e3a6SLisandro Dalcin } PC_HYPRE; 7116d9e3a6SLisandro Dalcin 72d2128fa2SBarry Smith #undef __FUNCT__ 73d2128fa2SBarry Smith #define __FUNCT__ "PCHYPREGetSolver" 74d2128fa2SBarry Smith PetscErrorCode PCHYPREGetSolver(PC pc,HYPRE_Solver *hsolver) 75d2128fa2SBarry Smith { 76d2128fa2SBarry Smith PC_HYPRE *jac = (PC_HYPRE*)pc->data; 77d2128fa2SBarry Smith 78d2128fa2SBarry Smith PetscFunctionBegin; 79d2128fa2SBarry Smith *hsolver = jac->hsolver; 80d2128fa2SBarry Smith PetscFunctionReturn(0); 81d2128fa2SBarry Smith } 8216d9e3a6SLisandro Dalcin 8316d9e3a6SLisandro Dalcin #undef __FUNCT__ 8416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetUp_HYPRE" 8516d9e3a6SLisandro Dalcin static PetscErrorCode PCSetUp_HYPRE(PC pc) 8616d9e3a6SLisandro Dalcin { 8716d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 8816d9e3a6SLisandro Dalcin PetscErrorCode ierr; 8916d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 9016d9e3a6SLisandro Dalcin HYPRE_ParVector bv,xv; 9116d9e3a6SLisandro Dalcin PetscInt bs; 9216d9e3a6SLisandro Dalcin 9316d9e3a6SLisandro Dalcin PetscFunctionBegin; 9416d9e3a6SLisandro Dalcin if (!jac->hypre_type) { 9502a17cd4SBarry Smith ierr = PCHYPRESetType(pc,"boomeramg");CHKERRQ(ierr); 9616d9e3a6SLisandro Dalcin } 975f5c5b43SBarry Smith 985f5c5b43SBarry Smith if (pc->setupcalled) { 995f5c5b43SBarry Smith /* always destroy the old matrix and create a new memory; 1005f5c5b43SBarry Smith hope this does not churn the memory too much. The problem 1015f5c5b43SBarry Smith is I do not know if it is possible to put the matrix back to 1025f5c5b43SBarry Smith its initial state so that we can directly copy the values 1035f5c5b43SBarry Smith the second time through. */ 104fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->ij)); 1055f5c5b43SBarry Smith jac->ij = 0; 10616d9e3a6SLisandro Dalcin } 1075f5c5b43SBarry Smith 10816d9e3a6SLisandro Dalcin if (!jac->ij) { /* create the matrix the first time through */ 10916d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCreate(pc->pmat,&jac->ij);CHKERRQ(ierr); 11016d9e3a6SLisandro Dalcin } 11116d9e3a6SLisandro Dalcin if (!jac->b) { /* create the vectors the first time through */ 11216d9e3a6SLisandro Dalcin Vec x,b; 1132a7a6963SBarry Smith ierr = MatCreateVecs(pc->pmat,&x,&b);CHKERRQ(ierr); 11416d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(x,&jac->x);CHKERRQ(ierr); 11516d9e3a6SLisandro Dalcin ierr = VecHYPRE_IJVectorCreate(b,&jac->b);CHKERRQ(ierr); 1166bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 1176bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 11816d9e3a6SLisandro Dalcin } 1195f5c5b43SBarry Smith 12016d9e3a6SLisandro Dalcin /* special case for BoomerAMG */ 12116d9e3a6SLisandro Dalcin if (jac->setup == HYPRE_BoomerAMGSetup) { 12216d9e3a6SLisandro Dalcin ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr); 1232fa5cd67SKarl Rupp if (bs > 1) PetscStackCallStandard(HYPRE_BoomerAMGSetNumFunctions,(jac->hsolver,bs)); 12416d9e3a6SLisandro Dalcin }; 12516d9e3a6SLisandro Dalcin ierr = MatHYPRE_IJMatrixCopy(pc->pmat,jac->ij);CHKERRQ(ierr); 126fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat)); 127fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&bv)); 128fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&xv)); 129fd3f9acdSBarry Smith PetscStackCall("HYPRE_SetupXXX",ierr = (*jac->setup)(jac->hsolver,hmat,bv,xv);CHKERRQ(ierr);); 13016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 13116d9e3a6SLisandro Dalcin } 13216d9e3a6SLisandro Dalcin 13316d9e3a6SLisandro Dalcin /* 13416d9e3a6SLisandro Dalcin Replaces the address where the HYPRE vector points to its data with the address of 13516d9e3a6SLisandro Dalcin PETSc's data. Saves the old address so it can be reset when we are finished with it. 13616d9e3a6SLisandro Dalcin Allows use to get the data into a HYPRE vector without the cost of memcopies 13716d9e3a6SLisandro Dalcin */ 13816d9e3a6SLisandro Dalcin #define HYPREReplacePointer(b,newvalue,savedvalue) { \ 13916d9e3a6SLisandro Dalcin hypre_ParVector *par_vector = (hypre_ParVector*)hypre_IJVectorObject(((hypre_IJVector*)b)); \ 14016d9e3a6SLisandro Dalcin hypre_Vector *local_vector = hypre_ParVectorLocalVector(par_vector); \ 14116d9e3a6SLisandro Dalcin savedvalue = local_vector->data; \ 1420ad7597dSKarl Rupp local_vector->data = newvalue; \ 1430ad7597dSKarl Rupp } 14416d9e3a6SLisandro Dalcin 14516d9e3a6SLisandro Dalcin #undef __FUNCT__ 14616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApply_HYPRE" 14716d9e3a6SLisandro Dalcin static PetscErrorCode PCApply_HYPRE(PC pc,Vec b,Vec x) 14816d9e3a6SLisandro Dalcin { 14916d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 15016d9e3a6SLisandro Dalcin PetscErrorCode ierr; 15116d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 15216d9e3a6SLisandro Dalcin PetscScalar *bv,*xv; 15316d9e3a6SLisandro Dalcin HYPRE_ParVector jbv,jxv; 15416d9e3a6SLisandro Dalcin PetscScalar *sbv,*sxv; 1554ddd07fcSJed Brown PetscInt hierr; 15616d9e3a6SLisandro Dalcin 15716d9e3a6SLisandro Dalcin PetscFunctionBegin; 158dff31646SBarry Smith ierr = PetscCitationsRegister(hypreCitation,&cite);CHKERRQ(ierr); 15916d9e3a6SLisandro Dalcin if (!jac->applyrichardson) {ierr = VecSet(x,0.0);CHKERRQ(ierr);} 16016d9e3a6SLisandro Dalcin ierr = VecGetArray(b,&bv);CHKERRQ(ierr); 16116d9e3a6SLisandro Dalcin ierr = VecGetArray(x,&xv);CHKERRQ(ierr); 16216d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,bv,sbv); 16316d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,xv,sxv); 16416d9e3a6SLisandro Dalcin 165fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat)); 166fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&jbv)); 167fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&jxv)); 168fd3f9acdSBarry Smith PetscStackCall("Hypre solve",hierr = (*jac->solve)(jac->hsolver,hmat,jbv,jxv); 16965e19b50SBarry Smith if (hierr && hierr != HYPRE_ERROR_CONV) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 170fd3f9acdSBarry Smith if (hierr) hypre__global_error = 0;); 17116d9e3a6SLisandro Dalcin 17216d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,sbv,bv); 17316d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,sxv,xv); 17416d9e3a6SLisandro Dalcin ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr); 17516d9e3a6SLisandro Dalcin ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr); 17616d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 17716d9e3a6SLisandro Dalcin } 17816d9e3a6SLisandro Dalcin 17916d9e3a6SLisandro Dalcin #undef __FUNCT__ 18016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCDestroy_HYPRE" 18116d9e3a6SLisandro Dalcin static PetscErrorCode PCDestroy_HYPRE(PC pc) 18216d9e3a6SLisandro Dalcin { 18316d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 18416d9e3a6SLisandro Dalcin PetscErrorCode ierr; 18516d9e3a6SLisandro Dalcin 18616d9e3a6SLisandro Dalcin PetscFunctionBegin; 187fd3f9acdSBarry Smith if (jac->ij) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(jac->ij)); 188fd3f9acdSBarry Smith if (jac->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->b)); 189fd3f9acdSBarry Smith if (jac->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(jac->x)); 190226b0620SJed Brown if (jac->destroy) PetscStackCall("HYPRE_DestroyXXX",ierr = (*jac->destroy)(jac->hsolver);CHKERRQ(ierr);); 191503cfb0cSBarry Smith ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr); 19216d9e3a6SLisandro Dalcin if (jac->comm_hypre != MPI_COMM_NULL) { ierr = MPI_Comm_free(&(jac->comm_hypre));CHKERRQ(ierr);} 193c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 19416d9e3a6SLisandro Dalcin 19516d9e3a6SLisandro Dalcin ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr); 196bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetType_C",NULL);CHKERRQ(ierr); 197bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHYPREGetType_C",NULL);CHKERRQ(ierr); 19816d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 19916d9e3a6SLisandro Dalcin } 20016d9e3a6SLisandro Dalcin 20116d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 20216d9e3a6SLisandro Dalcin #undef __FUNCT__ 20316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_Pilut" 2043931853cSBarry Smith static PetscErrorCode PCSetFromOptions_HYPRE_Pilut(PetscOptionsObjectType *PetscOptionsObject,PC pc) 20516d9e3a6SLisandro Dalcin { 20616d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 20716d9e3a6SLisandro Dalcin PetscErrorCode ierr; 208ace3abfcSBarry Smith PetscBool flag; 20916d9e3a6SLisandro Dalcin 21016d9e3a6SLisandro Dalcin PetscFunctionBegin; 211e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HYPRE Pilut Options");CHKERRQ(ierr); 21216d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_pilut_maxiter","Number of iterations","None",jac->maxiter,&jac->maxiter,&flag);CHKERRQ(ierr); 213fd3f9acdSBarry Smith if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetMaxIter,(jac->hsolver,jac->maxiter)); 21416d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_pilut_tol","Drop tolerance","None",jac->tol,&jac->tol,&flag);CHKERRQ(ierr); 215fd3f9acdSBarry Smith if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetDropTolerance,(jac->hsolver,jac->tol)); 21616d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_pilut_factorrowsize","FactorRowSize","None",jac->factorrowsize,&jac->factorrowsize,&flag);CHKERRQ(ierr); 217fd3f9acdSBarry Smith if (flag) PetscStackCallStandard(HYPRE_ParCSRPilutSetFactorRowSize,(jac->hsolver,jac->factorrowsize)); 21816d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 21916d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 22016d9e3a6SLisandro Dalcin } 22116d9e3a6SLisandro Dalcin 22216d9e3a6SLisandro Dalcin #undef __FUNCT__ 22316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_Pilut" 22416d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_Pilut(PC pc,PetscViewer viewer) 22516d9e3a6SLisandro Dalcin { 22616d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 22716d9e3a6SLisandro Dalcin PetscErrorCode ierr; 228ace3abfcSBarry Smith PetscBool iascii; 22916d9e3a6SLisandro Dalcin 23016d9e3a6SLisandro Dalcin PetscFunctionBegin; 231251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 23216d9e3a6SLisandro Dalcin if (iascii) { 23316d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut preconditioning\n");CHKERRQ(ierr); 23416d9e3a6SLisandro Dalcin if (jac->maxiter != PETSC_DEFAULT) { 23516d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: maximum number of iterations %d\n",jac->maxiter);CHKERRQ(ierr); 23616d9e3a6SLisandro Dalcin } else { 23716d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default maximum number of iterations \n");CHKERRQ(ierr); 23816d9e3a6SLisandro Dalcin } 23916d9e3a6SLisandro Dalcin if (jac->tol != PETSC_DEFAULT) { 24057622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: drop tolerance %g\n",(double)jac->tol);CHKERRQ(ierr); 24116d9e3a6SLisandro Dalcin } else { 24216d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default drop tolerance \n");CHKERRQ(ierr); 24316d9e3a6SLisandro Dalcin } 24416d9e3a6SLisandro Dalcin if (jac->factorrowsize != PETSC_DEFAULT) { 24516d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: factor row size %d\n",jac->factorrowsize);CHKERRQ(ierr); 24616d9e3a6SLisandro Dalcin } else { 24716d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE Pilut: default factor row size \n");CHKERRQ(ierr); 24816d9e3a6SLisandro Dalcin } 24916d9e3a6SLisandro Dalcin } 25016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 25116d9e3a6SLisandro Dalcin } 25216d9e3a6SLisandro Dalcin 25316d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 25416d9e3a6SLisandro Dalcin 25516d9e3a6SLisandro Dalcin #undef __FUNCT__ 25616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyTranspose_HYPRE_BoomerAMG" 25716d9e3a6SLisandro Dalcin static PetscErrorCode PCApplyTranspose_HYPRE_BoomerAMG(PC pc,Vec b,Vec x) 25816d9e3a6SLisandro Dalcin { 25916d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 26016d9e3a6SLisandro Dalcin PetscErrorCode ierr; 26116d9e3a6SLisandro Dalcin HYPRE_ParCSRMatrix hmat; 26216d9e3a6SLisandro Dalcin PetscScalar *bv,*xv; 26316d9e3a6SLisandro Dalcin HYPRE_ParVector jbv,jxv; 26416d9e3a6SLisandro Dalcin PetscScalar *sbv,*sxv; 2654ddd07fcSJed Brown PetscInt hierr; 26616d9e3a6SLisandro Dalcin 26716d9e3a6SLisandro Dalcin PetscFunctionBegin; 268dff31646SBarry Smith ierr = PetscCitationsRegister(hypreCitation,&cite);CHKERRQ(ierr); 26916d9e3a6SLisandro Dalcin ierr = VecSet(x,0.0);CHKERRQ(ierr); 27016d9e3a6SLisandro Dalcin ierr = VecGetArray(b,&bv);CHKERRQ(ierr); 27116d9e3a6SLisandro Dalcin ierr = VecGetArray(x,&xv);CHKERRQ(ierr); 27216d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,bv,sbv); 27316d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,xv,sxv); 27416d9e3a6SLisandro Dalcin 275fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJMatrixGetObject,(jac->ij,(void**)&hmat)); 276fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->b,(void**)&jbv)); 277fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_IJVectorGetObject,(jac->x,(void**)&jxv)); 27816d9e3a6SLisandro Dalcin 27916d9e3a6SLisandro Dalcin hierr = HYPRE_BoomerAMGSolveT(jac->hsolver,hmat,jbv,jxv); 28016d9e3a6SLisandro Dalcin /* error code of 1 in BoomerAMG merely means convergence not achieved */ 281e32f2f54SBarry Smith if (hierr && (hierr != 1)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HYPRE solver, error code %d",hierr); 28216d9e3a6SLisandro Dalcin if (hierr) hypre__global_error = 0; 28316d9e3a6SLisandro Dalcin 28416d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->b,sbv,bv); 28516d9e3a6SLisandro Dalcin HYPREReplacePointer(jac->x,sxv,xv); 28616d9e3a6SLisandro Dalcin ierr = VecRestoreArray(x,&xv);CHKERRQ(ierr); 28716d9e3a6SLisandro Dalcin ierr = VecRestoreArray(b,&bv);CHKERRQ(ierr); 28816d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 28916d9e3a6SLisandro Dalcin } 29016d9e3a6SLisandro Dalcin 291a669f990SJed Brown /* static array length */ 292a669f990SJed Brown #define ALEN(a) (sizeof(a)/sizeof((a)[0])) 293a669f990SJed Brown 29416d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGCycleType[] = {"","V","W"}; 2950f1074feSSatish Balay static const char *HYPREBoomerAMGCoarsenType[] = {"CLJP","Ruge-Stueben","","modifiedRuge-Stueben","","","Falgout", "", "PMIS", "", "HMIS"}; 29616d9e3a6SLisandro Dalcin static const char *HYPREBoomerAMGMeasureType[] = {"local","global"}; 29765de4495SJed Brown /* The following corresponds to HYPRE_BoomerAMGSetRelaxType which has many missing numbers in the enum */ 29865de4495SJed Brown static const char *HYPREBoomerAMGRelaxType[] = {"Jacobi","sequential-Gauss-Seidel","seqboundary-Gauss-Seidel","SOR/Jacobi","backward-SOR/Jacobi", 29965de4495SJed Brown "" /* [5] hybrid chaotic Gauss-Seidel (works only with OpenMP) */,"symmetric-SOR/Jacobi", 30065de4495SJed Brown "" /* 7 */,"l1scaled-SOR/Jacobi","Gaussian-elimination", 30165de4495SJed Brown "" /* 10 */, "" /* 11 */, "" /* 12 */, "" /* 13 */, "" /* 14 */, 30265de4495SJed Brown "CG" /* non-stationary */,"Chebyshev","FCF-Jacobi","l1scaled-Jacobi"}; 3030f1074feSSatish Balay static const char *HYPREBoomerAMGInterpType[] = {"classical", "", "", "direct", "multipass", "multipass-wts", "ext+i", 3040f1074feSSatish Balay "ext+i-cc", "standard", "standard-wts", "", "", "FF", "FF1"}; 30516d9e3a6SLisandro Dalcin #undef __FUNCT__ 30616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_BoomerAMG" 3073931853cSBarry Smith static PetscErrorCode PCSetFromOptions_HYPRE_BoomerAMG(PetscOptionsObjectType *PetscOptionsObject,PC pc) 30816d9e3a6SLisandro Dalcin { 30916d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 31016d9e3a6SLisandro Dalcin PetscErrorCode ierr; 3114ddd07fcSJed Brown PetscInt n,indx,level; 312ace3abfcSBarry Smith PetscBool flg, tmp_truth; 31316d9e3a6SLisandro Dalcin double tmpdbl, twodbl[2]; 31416d9e3a6SLisandro Dalcin 31516d9e3a6SLisandro Dalcin PetscFunctionBegin; 316e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HYPRE BoomerAMG Options");CHKERRQ(ierr); 3174336a9eeSBarry Smith ierr = PetscOptionsEList("-pc_hypre_boomeramg_cycle_type","Cycle type","None",HYPREBoomerAMGCycleType+1,2,HYPREBoomerAMGCycleType[jac->cycletype],&indx,&flg);CHKERRQ(ierr); 31816d9e3a6SLisandro Dalcin if (flg) { 3194336a9eeSBarry Smith jac->cycletype = indx+1; 320fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype)); 32116d9e3a6SLisandro Dalcin } 32216d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_levels","Number of levels (of grids) allowed","None",jac->maxlevels,&jac->maxlevels,&flg);CHKERRQ(ierr); 32316d9e3a6SLisandro Dalcin if (flg) { 324ce94432eSBarry Smith if (jac->maxlevels < 2) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %d must be at least two",jac->maxlevels); 325fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels)); 32616d9e3a6SLisandro Dalcin } 32716d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_max_iter","Maximum iterations used PER hypre call","None",jac->maxiter,&jac->maxiter,&flg);CHKERRQ(ierr); 32816d9e3a6SLisandro Dalcin if (flg) { 329ce94432eSBarry Smith if (jac->maxiter < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of iterations %d must be at least one",jac->maxiter); 330fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter)); 33116d9e3a6SLisandro Dalcin } 3320f1074feSSatish 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); 33316d9e3a6SLisandro Dalcin if (flg) { 33457622a8eSBarry Smith if (jac->tol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Tolerance %g must be greater than or equal to zero",(double)jac->tol); 335fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol)); 33616d9e3a6SLisandro Dalcin } 33716d9e3a6SLisandro Dalcin 3380f1074feSSatish Balay ierr = PetscOptionsScalar("-pc_hypre_boomeramg_truncfactor","Truncation factor for interpolation (0=no truncation)","None",jac->truncfactor,&jac->truncfactor,&flg);CHKERRQ(ierr); 33916d9e3a6SLisandro Dalcin if (flg) { 34057622a8eSBarry Smith if (jac->truncfactor < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Truncation factor %g must be great than or equal zero",(double)jac->truncfactor); 341fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor)); 34216d9e3a6SLisandro Dalcin } 34316d9e3a6SLisandro Dalcin 3440f1074feSSatish 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); 3450f1074feSSatish Balay if (flg) { 34657622a8eSBarry Smith if (jac->pmax < 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"P_max %g must be greater than or equal to zero",(double)jac->pmax); 347fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax)); 3480f1074feSSatish Balay } 3490f1074feSSatish Balay 3500f1074feSSatish Balay ierr = PetscOptionsInt("-pc_hypre_boomeramg_agg_nl","Number of levels of aggressive coarsening","None",jac->agg_nl,&jac->agg_nl,&flg);CHKERRQ(ierr); 3510f1074feSSatish Balay if (flg) { 35257622a8eSBarry Smith if (jac->agg_nl < 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of levels %g must be greater than or equal to zero",(double)jac->agg_nl); 3530f1074feSSatish Balay 354fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl)); 3550f1074feSSatish Balay } 3560f1074feSSatish Balay 3570f1074feSSatish Balay 3580f1074feSSatish 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); 3590f1074feSSatish Balay if (flg) { 36057622a8eSBarry Smith if (jac->agg_num_paths < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of paths %g must be greater than or equal to 1",(double)jac->agg_num_paths); 3610f1074feSSatish Balay 362fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths)); 3630f1074feSSatish Balay } 3640f1074feSSatish Balay 3650f1074feSSatish Balay 36616d9e3a6SLisandro Dalcin ierr = PetscOptionsScalar("-pc_hypre_boomeramg_strong_threshold","Threshold for being strongly connected","None",jac->strongthreshold,&jac->strongthreshold,&flg);CHKERRQ(ierr); 36716d9e3a6SLisandro Dalcin if (flg) { 36857622a8eSBarry Smith if (jac->strongthreshold < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Strong threshold %g must be great than or equal zero",(double)jac->strongthreshold); 369fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold)); 37016d9e3a6SLisandro Dalcin } 37116d9e3a6SLisandro Dalcin ierr = PetscOptionsScalar("-pc_hypre_boomeramg_max_row_sum","Maximum row sum","None",jac->maxrowsum,&jac->maxrowsum,&flg);CHKERRQ(ierr); 37216d9e3a6SLisandro Dalcin if (flg) { 37357622a8eSBarry Smith if (jac->maxrowsum < 0.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %g must be greater than zero",(double)jac->maxrowsum); 37457622a8eSBarry Smith if (jac->maxrowsum > 1.0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Maximum row sum %g must be less than or equal one",(double)jac->maxrowsum); 375fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum)); 37616d9e3a6SLisandro Dalcin } 37716d9e3a6SLisandro Dalcin 37816d9e3a6SLisandro Dalcin /* Grid sweeps */ 3790f1074feSSatish 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); 38016d9e3a6SLisandro Dalcin if (flg) { 381fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver,indx)); 38216d9e3a6SLisandro Dalcin /* modify the jac structure so we can view the updated options with PC_View */ 38316d9e3a6SLisandro Dalcin jac->gridsweeps[0] = indx; 3840f1074feSSatish Balay jac->gridsweeps[1] = indx; 3850f1074feSSatish Balay /*defaults coarse to 1 */ 3860f1074feSSatish Balay jac->gridsweeps[2] = 1; 38716d9e3a6SLisandro Dalcin } 3880f1074feSSatish Balay 3890f1074feSSatish Balay ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_down","Number of sweeps for the down cycles","None",jac->gridsweeps[0], &indx,&flg);CHKERRQ(ierr); 39016d9e3a6SLisandro Dalcin if (flg) { 391fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 1)); 3920f1074feSSatish Balay jac->gridsweeps[0] = indx; 39316d9e3a6SLisandro Dalcin } 39416d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_up","Number of sweeps for the up cycles","None",jac->gridsweeps[1],&indx,&flg);CHKERRQ(ierr); 39516d9e3a6SLisandro Dalcin if (flg) { 396fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 2)); 3970f1074feSSatish Balay jac->gridsweeps[1] = indx; 39816d9e3a6SLisandro Dalcin } 3990f1074feSSatish Balay ierr = PetscOptionsInt("-pc_hypre_boomeramg_grid_sweeps_coarse","Number of sweeps for the coarse level","None",jac->gridsweeps[2],&indx,&flg);CHKERRQ(ierr); 40016d9e3a6SLisandro Dalcin if (flg) { 401fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetCycleNumSweeps,(jac->hsolver,indx, 3)); 4020f1074feSSatish Balay jac->gridsweeps[2] = indx; 40316d9e3a6SLisandro Dalcin } 40416d9e3a6SLisandro Dalcin 40516d9e3a6SLisandro Dalcin /* Relax type */ 406a669f990SJed Brown ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_all","Relax type for the up and down cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr); 40716d9e3a6SLisandro Dalcin if (flg) { 4080f1074feSSatish Balay jac->relaxtype[0] = jac->relaxtype[1] = indx; 409fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, indx)); 4100f1074feSSatish Balay /* by default, coarse type set to 9 */ 4110f1074feSSatish Balay jac->relaxtype[2] = 9; 4120f1074feSSatish Balay 41316d9e3a6SLisandro Dalcin } 414a669f990SJed Brown ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_down","Relax type for the down cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr); 41516d9e3a6SLisandro Dalcin if (flg) { 41616d9e3a6SLisandro Dalcin jac->relaxtype[0] = indx; 417fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 1)); 41816d9e3a6SLisandro Dalcin } 419a669f990SJed Brown ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_up","Relax type for the up cycles","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[6],&indx,&flg);CHKERRQ(ierr); 42016d9e3a6SLisandro Dalcin if (flg) { 4210f1074feSSatish Balay jac->relaxtype[1] = indx; 422fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 2)); 42316d9e3a6SLisandro Dalcin } 424a669f990SJed Brown ierr = PetscOptionsEList("-pc_hypre_boomeramg_relax_type_coarse","Relax type on coarse grid","None",HYPREBoomerAMGRelaxType,ALEN(HYPREBoomerAMGRelaxType),HYPREBoomerAMGRelaxType[9],&indx,&flg);CHKERRQ(ierr); 42516d9e3a6SLisandro Dalcin if (flg) { 4260f1074feSSatish Balay jac->relaxtype[2] = indx; 427fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetCycleRelaxType,(jac->hsolver, indx, 3)); 42816d9e3a6SLisandro Dalcin } 42916d9e3a6SLisandro Dalcin 43016d9e3a6SLisandro Dalcin /* Relaxation Weight */ 43116d9e3a6SLisandro 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); 43216d9e3a6SLisandro Dalcin if (flg) { 433fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxWt,(jac->hsolver,tmpdbl)); 43416d9e3a6SLisandro Dalcin jac->relaxweight = tmpdbl; 43516d9e3a6SLisandro Dalcin } 43616d9e3a6SLisandro Dalcin 43716d9e3a6SLisandro Dalcin n = 2; 43816d9e3a6SLisandro Dalcin twodbl[0] = twodbl[1] = 1.0; 43916d9e3a6SLisandro 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); 44016d9e3a6SLisandro Dalcin if (flg) { 44116d9e3a6SLisandro Dalcin if (n == 2) { 44216d9e3a6SLisandro Dalcin indx = (int)PetscAbsReal(twodbl[1]); 443fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetLevelRelaxWt,(jac->hsolver,twodbl[0],indx)); 444ce94432eSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relax weight level: you must provide 2 values separated by a comma (and no space), you provided %d",n); 44516d9e3a6SLisandro Dalcin } 44616d9e3a6SLisandro Dalcin 44716d9e3a6SLisandro Dalcin /* Outer relaxation Weight */ 44816d9e3a6SLisandro 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); 44916d9e3a6SLisandro Dalcin if (flg) { 450fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetOuterWt,(jac->hsolver, tmpdbl)); 45116d9e3a6SLisandro Dalcin jac->outerrelaxweight = tmpdbl; 45216d9e3a6SLisandro Dalcin } 45316d9e3a6SLisandro Dalcin 45416d9e3a6SLisandro Dalcin n = 2; 45516d9e3a6SLisandro Dalcin twodbl[0] = twodbl[1] = 1.0; 45616d9e3a6SLisandro 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); 45716d9e3a6SLisandro Dalcin if (flg) { 45816d9e3a6SLisandro Dalcin if (n == 2) { 45916d9e3a6SLisandro Dalcin indx = (int)PetscAbsReal(twodbl[1]); 460fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetLevelOuterWt,(jac->hsolver, twodbl[0], indx)); 461ce94432eSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relax weight outer level: You must provide 2 values separated by a comma (and no space), you provided %d",n); 46216d9e3a6SLisandro Dalcin } 46316d9e3a6SLisandro Dalcin 46416d9e3a6SLisandro Dalcin /* the Relax Order */ 465acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_hypre_boomeramg_no_CF", "Do not use CF-relaxation", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr); 46616d9e3a6SLisandro Dalcin 4678afaa268SBarry Smith if (flg && tmp_truth) { 46816d9e3a6SLisandro Dalcin jac->relaxorder = 0; 469fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder)); 47016d9e3a6SLisandro Dalcin } 471a669f990SJed Brown ierr = PetscOptionsEList("-pc_hypre_boomeramg_measure_type","Measure type","None",HYPREBoomerAMGMeasureType,ALEN(HYPREBoomerAMGMeasureType),HYPREBoomerAMGMeasureType[0],&indx,&flg);CHKERRQ(ierr); 47216d9e3a6SLisandro Dalcin if (flg) { 47316d9e3a6SLisandro Dalcin jac->measuretype = indx; 474fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype)); 47516d9e3a6SLisandro Dalcin } 4760f1074feSSatish Balay /* update list length 3/07 */ 477a669f990SJed Brown ierr = PetscOptionsEList("-pc_hypre_boomeramg_coarsen_type","Coarsen type","None",HYPREBoomerAMGCoarsenType,ALEN(HYPREBoomerAMGCoarsenType),HYPREBoomerAMGCoarsenType[6],&indx,&flg);CHKERRQ(ierr); 47816d9e3a6SLisandro Dalcin if (flg) { 47916d9e3a6SLisandro Dalcin jac->coarsentype = indx; 480fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype)); 48116d9e3a6SLisandro Dalcin } 4820f1074feSSatish Balay 4830f1074feSSatish Balay /* new 3/07 */ 484a669f990SJed Brown ierr = PetscOptionsEList("-pc_hypre_boomeramg_interp_type","Interpolation type","None",HYPREBoomerAMGInterpType,ALEN(HYPREBoomerAMGInterpType),HYPREBoomerAMGInterpType[0],&indx,&flg);CHKERRQ(ierr); 4850f1074feSSatish Balay if (flg) { 4860f1074feSSatish Balay jac->interptype = indx; 487fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype)); 4880f1074feSSatish Balay } 4890f1074feSSatish Balay 490b96a4a96SBarry Smith ierr = PetscOptionsName("-pc_hypre_boomeramg_print_statistics","Print statistics","None",&flg);CHKERRQ(ierr); 49116d9e3a6SLisandro Dalcin if (flg) { 492b96a4a96SBarry Smith level = 3; 4930298fd71SBarry Smith ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_statistics","Print statistics","None",level,&level,NULL);CHKERRQ(ierr); 4942fa5cd67SKarl Rupp 495b96a4a96SBarry Smith jac->printstatistics = PETSC_TRUE; 496fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetPrintLevel,(jac->hsolver,level)); 4972ae77aedSBarry Smith } 4982ae77aedSBarry Smith 499b96a4a96SBarry Smith ierr = PetscOptionsName("-pc_hypre_boomeramg_print_debug","Print debug information","None",&flg);CHKERRQ(ierr); 5002ae77aedSBarry Smith if (flg) { 501b96a4a96SBarry Smith level = 3; 5020298fd71SBarry Smith ierr = PetscOptionsInt("-pc_hypre_boomeramg_print_debug","Print debug information","None",level,&level,NULL);CHKERRQ(ierr); 5032fa5cd67SKarl Rupp 504b96a4a96SBarry Smith jac->printstatistics = PETSC_TRUE; 505fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetDebugFlag,(jac->hsolver,level)); 50616d9e3a6SLisandro Dalcin } 5078f87f92bSBarry Smith 5088afaa268SBarry Smith ierr = PetscOptionsBool("-pc_hypre_boomeramg_nodal_coarsen", "HYPRE_BoomerAMGSetNodal()", "None", jac->nodal_coarsen ? PETSC_TRUE : PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr); 5098f87f92bSBarry Smith if (flg && tmp_truth) { 5108f87f92bSBarry Smith jac->nodal_coarsen = 1; 511fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetNodal,(jac->hsolver,1)); 5128f87f92bSBarry Smith } 5138f87f92bSBarry Smith 514acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None", PETSC_FALSE, &tmp_truth, &flg);CHKERRQ(ierr); 5158f87f92bSBarry Smith if (flg && tmp_truth) { 5168f87f92bSBarry Smith PetscInt tmp_int; 5178f87f92bSBarry Smith ierr = PetscOptionsInt("-pc_hypre_boomeramg_nodal_relaxation", "Nodal relaxation via Schwarz", "None",jac->nodal_relax_levels,&tmp_int,&flg);CHKERRQ(ierr); 5188f87f92bSBarry Smith if (flg) jac->nodal_relax_levels = tmp_int; 519fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothType,(jac->hsolver,6)); 520fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetDomainType,(jac->hsolver,1)); 521fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetOverlap,(jac->hsolver,0)); 522fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetSmoothNumLevels,(jac->hsolver,jac->nodal_relax_levels)); 5238f87f92bSBarry Smith } 5248f87f92bSBarry Smith 52516d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 52616d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 52716d9e3a6SLisandro Dalcin } 52816d9e3a6SLisandro Dalcin 52916d9e3a6SLisandro Dalcin #undef __FUNCT__ 53016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCApplyRichardson_HYPRE_BoomerAMG" 531ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_HYPRE_BoomerAMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason) 53216d9e3a6SLisandro Dalcin { 53316d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 53416d9e3a6SLisandro Dalcin PetscErrorCode ierr; 5354ddd07fcSJed Brown PetscInt oits; 53616d9e3a6SLisandro Dalcin 53716d9e3a6SLisandro Dalcin PetscFunctionBegin; 538dff31646SBarry Smith ierr = PetscCitationsRegister(hypreCitation,&cite);CHKERRQ(ierr); 539fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,its*jac->maxiter)); 540fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,rtol)); 54116d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_TRUE; 54216d9e3a6SLisandro Dalcin ierr = PCApply_HYPRE(pc,b,y);CHKERRQ(ierr); 54316d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_FALSE; 5448b1f7689SBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGGetNumIterations,(jac->hsolver,(HYPRE_Int *)&oits)); 5454d0a8057SBarry Smith *outits = oits; 5464d0a8057SBarry Smith if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS; 5474d0a8057SBarry Smith else *reason = PCRICHARDSON_CONVERGED_RTOL; 548fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol)); 549fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter)); 55016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 55116d9e3a6SLisandro Dalcin } 55216d9e3a6SLisandro Dalcin 55316d9e3a6SLisandro Dalcin 55416d9e3a6SLisandro Dalcin #undef __FUNCT__ 55516d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_BoomerAMG" 55616d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_BoomerAMG(PC pc,PetscViewer viewer) 55716d9e3a6SLisandro Dalcin { 55816d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 55916d9e3a6SLisandro Dalcin PetscErrorCode ierr; 560ace3abfcSBarry Smith PetscBool iascii; 56116d9e3a6SLisandro Dalcin 56216d9e3a6SLisandro Dalcin PetscFunctionBegin; 563251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 56416d9e3a6SLisandro Dalcin if (iascii) { 56516d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG preconditioning\n");CHKERRQ(ierr); 56616d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Cycle type %s\n",HYPREBoomerAMGCycleType[jac->cycletype]);CHKERRQ(ierr); 56716d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum number of levels %d\n",jac->maxlevels);CHKERRQ(ierr); 56816d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum number of iterations PER hypre call %d\n",jac->maxiter);CHKERRQ(ierr); 56957622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Convergence tolerance PER hypre call %g\n",(double)jac->tol);CHKERRQ(ierr); 57057622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Threshold for strong coupling %g\n",(double)jac->strongthreshold);CHKERRQ(ierr); 57157622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation truncation factor %g\n",(double)jac->truncfactor);CHKERRQ(ierr); 5720f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation: max elements per row %d\n",jac->pmax);CHKERRQ(ierr); 5730f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Number of levels of aggressive coarsening %d\n",jac->agg_nl);CHKERRQ(ierr); 5740f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Number of paths for aggressive coarsening %d\n",jac->agg_num_paths);CHKERRQ(ierr); 5750f1074feSSatish Balay 57657622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Maximum row sums %g\n",(double)jac->maxrowsum);CHKERRQ(ierr); 57716d9e3a6SLisandro Dalcin 5780f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps down %d\n",jac->gridsweeps[0]);CHKERRQ(ierr); 5790f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps up %d\n",jac->gridsweeps[1]);CHKERRQ(ierr); 5800f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Sweeps on coarse %d\n",jac->gridsweeps[2]);CHKERRQ(ierr); 58116d9e3a6SLisandro Dalcin 5820f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax down %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[0]]);CHKERRQ(ierr); 5830f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax up %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[1]]);CHKERRQ(ierr); 5840f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax on coarse %s\n",HYPREBoomerAMGRelaxType[jac->relaxtype[2]]);CHKERRQ(ierr); 58516d9e3a6SLisandro Dalcin 58657622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Relax weight (all) %g\n",(double)jac->relaxweight);CHKERRQ(ierr); 58757622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Outer relax weight (all) %g\n",(double)jac->outerrelaxweight);CHKERRQ(ierr); 58816d9e3a6SLisandro Dalcin 58916d9e3a6SLisandro Dalcin if (jac->relaxorder) { 59016d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using CF-relaxation\n");CHKERRQ(ierr); 59116d9e3a6SLisandro Dalcin } else { 59216d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Not using CF-relaxation\n");CHKERRQ(ierr); 59316d9e3a6SLisandro Dalcin } 59416d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Measure type %s\n",HYPREBoomerAMGMeasureType[jac->measuretype]);CHKERRQ(ierr); 59516d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Coarsen type %s\n",HYPREBoomerAMGCoarsenType[jac->coarsentype]);CHKERRQ(ierr); 5960f1074feSSatish Balay ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Interpolation type %s\n",HYPREBoomerAMGInterpType[jac->interptype]);CHKERRQ(ierr); 5978f87f92bSBarry Smith if (jac->nodal_coarsen) { 5988f87f92bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal coarsening (with HYPRE_BOOMERAMGSetNodal())\n");CHKERRQ(ierr); 5998f87f92bSBarry Smith } 6008f87f92bSBarry Smith if (jac->nodal_relax) { 6018f87f92bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE BoomerAMG: Using nodal relaxation via Schwarz smoothing on levels %d\n",jac->nodal_relax_levels);CHKERRQ(ierr); 6028f87f92bSBarry Smith } 60316d9e3a6SLisandro Dalcin } 60416d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 60516d9e3a6SLisandro Dalcin } 60616d9e3a6SLisandro Dalcin 60716d9e3a6SLisandro Dalcin /* --------------------------------------------------------------------------------------------*/ 60816d9e3a6SLisandro Dalcin #undef __FUNCT__ 60916d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE_ParaSails" 6103931853cSBarry Smith static PetscErrorCode PCSetFromOptions_HYPRE_ParaSails(PetscOptionsObjectType *PetscOptionsObject,PC pc) 61116d9e3a6SLisandro Dalcin { 61216d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 61316d9e3a6SLisandro Dalcin PetscErrorCode ierr; 6144ddd07fcSJed Brown PetscInt indx; 615ace3abfcSBarry Smith PetscBool flag; 61616d9e3a6SLisandro Dalcin const char *symtlist[] = {"nonsymmetric","SPD","nonsymmetric,SPD"}; 61716d9e3a6SLisandro Dalcin 61816d9e3a6SLisandro Dalcin PetscFunctionBegin; 619e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HYPRE ParaSails Options");CHKERRQ(ierr); 62016d9e3a6SLisandro Dalcin ierr = PetscOptionsInt("-pc_hypre_parasails_nlevels","Number of number of levels","None",jac->nlevels,&jac->nlevels,0);CHKERRQ(ierr); 62116d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_thresh","Threshold","None",jac->threshhold,&jac->threshhold,&flag);CHKERRQ(ierr); 6222fa5cd67SKarl Rupp if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels)); 62316d9e3a6SLisandro Dalcin 62416d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_filter","filter","None",jac->filter,&jac->filter,&flag);CHKERRQ(ierr); 6252fa5cd67SKarl Rupp if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter)); 62616d9e3a6SLisandro Dalcin 62716d9e3a6SLisandro Dalcin ierr = PetscOptionsReal("-pc_hypre_parasails_loadbal","Load balance","None",jac->loadbal,&jac->loadbal,&flag);CHKERRQ(ierr); 6282fa5cd67SKarl Rupp if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal)); 62916d9e3a6SLisandro Dalcin 630acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_hypre_parasails_logging","Print info to screen","None",(PetscBool)jac->logging,(PetscBool*)&jac->logging,&flag);CHKERRQ(ierr); 6312fa5cd67SKarl Rupp if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging)); 63216d9e3a6SLisandro Dalcin 633acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_hypre_parasails_reuse","Reuse nonzero pattern in preconditioner","None",(PetscBool)jac->ruse,(PetscBool*)&jac->ruse,&flag);CHKERRQ(ierr); 6342fa5cd67SKarl Rupp if (flag) PetscStackCallStandard(HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse)); 63516d9e3a6SLisandro Dalcin 636a669f990SJed Brown ierr = PetscOptionsEList("-pc_hypre_parasails_sym","Symmetry of matrix and preconditioner","None",symtlist,ALEN(symtlist),symtlist[0],&indx,&flag);CHKERRQ(ierr); 63716d9e3a6SLisandro Dalcin if (flag) { 63816d9e3a6SLisandro Dalcin jac->symt = indx; 639fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt)); 64016d9e3a6SLisandro Dalcin } 64116d9e3a6SLisandro Dalcin 64216d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 64316d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 64416d9e3a6SLisandro Dalcin } 64516d9e3a6SLisandro Dalcin 64616d9e3a6SLisandro Dalcin #undef __FUNCT__ 64716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCView_HYPRE_ParaSails" 64816d9e3a6SLisandro Dalcin static PetscErrorCode PCView_HYPRE_ParaSails(PC pc,PetscViewer viewer) 64916d9e3a6SLisandro Dalcin { 65016d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 65116d9e3a6SLisandro Dalcin PetscErrorCode ierr; 652ace3abfcSBarry Smith PetscBool iascii; 65316d9e3a6SLisandro Dalcin const char *symt = 0;; 65416d9e3a6SLisandro Dalcin 65516d9e3a6SLisandro Dalcin PetscFunctionBegin; 656251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 65716d9e3a6SLisandro Dalcin if (iascii) { 65816d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails preconditioning\n");CHKERRQ(ierr); 65916d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: nlevels %d\n",jac->nlevels);CHKERRQ(ierr); 66057622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: threshold %g\n",(double)jac->threshhold);CHKERRQ(ierr); 66157622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: filter %g\n",(double)jac->filter);CHKERRQ(ierr); 66257622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: load balance %g\n",(double)jac->loadbal);CHKERRQ(ierr); 663ace3abfcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: reuse nonzero structure %s\n",PetscBools[jac->ruse]);CHKERRQ(ierr); 664ace3abfcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: print info to screen %s\n",PetscBools[jac->logging]);CHKERRQ(ierr); 6652fa5cd67SKarl Rupp if (!jac->symt) symt = "nonsymmetric matrix and preconditioner"; 6662fa5cd67SKarl Rupp else if (jac->symt == 1) symt = "SPD matrix and preconditioner"; 6672fa5cd67SKarl Rupp else if (jac->symt == 2) symt = "nonsymmetric matrix but SPD preconditioner"; 668ce94432eSBarry Smith else SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Unknown HYPRE ParaSails symmetric option %d",jac->symt); 66916d9e3a6SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer," HYPRE ParaSails: %s\n",symt);CHKERRQ(ierr); 67016d9e3a6SLisandro Dalcin } 67116d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 67216d9e3a6SLisandro Dalcin } 67316d9e3a6SLisandro Dalcin /* ---------------------------------------------------------------------------------*/ 67416d9e3a6SLisandro Dalcin 67516d9e3a6SLisandro Dalcin #undef __FUNCT__ 67616d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType_HYPRE" 677f7a08781SBarry Smith static PetscErrorCode PCHYPREGetType_HYPRE(PC pc,const char *name[]) 67816d9e3a6SLisandro Dalcin { 67916d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 68016d9e3a6SLisandro Dalcin 68116d9e3a6SLisandro Dalcin PetscFunctionBegin; 68216d9e3a6SLisandro Dalcin *name = jac->hypre_type; 68316d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 68416d9e3a6SLisandro Dalcin } 68516d9e3a6SLisandro Dalcin 68616d9e3a6SLisandro Dalcin #undef __FUNCT__ 68716d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType_HYPRE" 688f7a08781SBarry Smith static PetscErrorCode PCHYPRESetType_HYPRE(PC pc,const char name[]) 68916d9e3a6SLisandro Dalcin { 69016d9e3a6SLisandro Dalcin PC_HYPRE *jac = (PC_HYPRE*)pc->data; 69116d9e3a6SLisandro Dalcin PetscErrorCode ierr; 692ace3abfcSBarry Smith PetscBool flag; 69316d9e3a6SLisandro Dalcin 69416d9e3a6SLisandro Dalcin PetscFunctionBegin; 69516d9e3a6SLisandro Dalcin if (jac->hypre_type) { 69616d9e3a6SLisandro Dalcin ierr = PetscStrcmp(jac->hypre_type,name,&flag);CHKERRQ(ierr); 697ce94432eSBarry Smith if (!flag) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Cannot reset the HYPRE preconditioner type once it has been set"); 69816d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 69916d9e3a6SLisandro Dalcin } else { 70016d9e3a6SLisandro Dalcin ierr = PetscStrallocpy(name, &jac->hypre_type);CHKERRQ(ierr); 70116d9e3a6SLisandro Dalcin } 70216d9e3a6SLisandro Dalcin 70316d9e3a6SLisandro Dalcin jac->maxiter = PETSC_DEFAULT; 70416d9e3a6SLisandro Dalcin jac->tol = PETSC_DEFAULT; 70516d9e3a6SLisandro Dalcin jac->printstatistics = PetscLogPrintInfo; 70616d9e3a6SLisandro Dalcin 70716d9e3a6SLisandro Dalcin ierr = PetscStrcmp("pilut",jac->hypre_type,&flag);CHKERRQ(ierr); 70816d9e3a6SLisandro Dalcin if (flag) { 709fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_ParCSRPilutCreate,(jac->comm_hypre,&jac->hsolver)); 71016d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_Pilut; 71116d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_Pilut; 71216d9e3a6SLisandro Dalcin jac->destroy = HYPRE_ParCSRPilutDestroy; 71316d9e3a6SLisandro Dalcin jac->setup = HYPRE_ParCSRPilutSetup; 71416d9e3a6SLisandro Dalcin jac->solve = HYPRE_ParCSRPilutSolve; 71516d9e3a6SLisandro Dalcin jac->factorrowsize = PETSC_DEFAULT; 71616d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 71716d9e3a6SLisandro Dalcin } 71816d9e3a6SLisandro Dalcin ierr = PetscStrcmp("parasails",jac->hypre_type,&flag);CHKERRQ(ierr); 71916d9e3a6SLisandro Dalcin if (flag) { 720fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_ParaSailsCreate,(jac->comm_hypre,&jac->hsolver)); 72116d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_ParaSails; 72216d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_ParaSails; 72316d9e3a6SLisandro Dalcin jac->destroy = HYPRE_ParaSailsDestroy; 72416d9e3a6SLisandro Dalcin jac->setup = HYPRE_ParaSailsSetup; 72516d9e3a6SLisandro Dalcin jac->solve = HYPRE_ParaSailsSolve; 72616d9e3a6SLisandro Dalcin /* initialize */ 72716d9e3a6SLisandro Dalcin jac->nlevels = 1; 72816d9e3a6SLisandro Dalcin jac->threshhold = .1; 72916d9e3a6SLisandro Dalcin jac->filter = .1; 73016d9e3a6SLisandro Dalcin jac->loadbal = 0; 7312fa5cd67SKarl Rupp if (PetscLogPrintInfo) jac->logging = (int) PETSC_TRUE; 7322fa5cd67SKarl Rupp else jac->logging = (int) PETSC_FALSE; 7332fa5cd67SKarl Rupp 73416d9e3a6SLisandro Dalcin jac->ruse = (int) PETSC_FALSE; 73516d9e3a6SLisandro Dalcin jac->symt = 0; 736fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_ParaSailsSetParams,(jac->hsolver,jac->threshhold,jac->nlevels)); 737fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_ParaSailsSetFilter,(jac->hsolver,jac->filter)); 738fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_ParaSailsSetLoadbal,(jac->hsolver,jac->loadbal)); 739fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_ParaSailsSetLogging,(jac->hsolver,jac->logging)); 740fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_ParaSailsSetReuse,(jac->hsolver,jac->ruse)); 741fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_ParaSailsSetSym,(jac->hsolver,jac->symt)); 74216d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 74316d9e3a6SLisandro Dalcin } 74416d9e3a6SLisandro Dalcin ierr = PetscStrcmp("boomeramg",jac->hypre_type,&flag);CHKERRQ(ierr); 74516d9e3a6SLisandro Dalcin if (flag) { 74616d9e3a6SLisandro Dalcin ierr = HYPRE_BoomerAMGCreate(&jac->hsolver); 74716d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE_BoomerAMG; 74816d9e3a6SLisandro Dalcin pc->ops->view = PCView_HYPRE_BoomerAMG; 74916d9e3a6SLisandro Dalcin pc->ops->applytranspose = PCApplyTranspose_HYPRE_BoomerAMG; 75016d9e3a6SLisandro Dalcin pc->ops->applyrichardson = PCApplyRichardson_HYPRE_BoomerAMG; 75116d9e3a6SLisandro Dalcin jac->destroy = HYPRE_BoomerAMGDestroy; 75216d9e3a6SLisandro Dalcin jac->setup = HYPRE_BoomerAMGSetup; 75316d9e3a6SLisandro Dalcin jac->solve = HYPRE_BoomerAMGSolve; 75416d9e3a6SLisandro Dalcin jac->applyrichardson = PETSC_FALSE; 75516d9e3a6SLisandro Dalcin /* these defaults match the hypre defaults */ 75616d9e3a6SLisandro Dalcin jac->cycletype = 1; 75716d9e3a6SLisandro Dalcin jac->maxlevels = 25; 75816d9e3a6SLisandro Dalcin jac->maxiter = 1; 7598f87f92bSBarry Smith jac->tol = 0.0; /* tolerance of zero indicates use as preconditioner (suppresses convergence errors) */ 76016d9e3a6SLisandro Dalcin jac->truncfactor = 0.0; 76116d9e3a6SLisandro Dalcin jac->strongthreshold = .25; 76216d9e3a6SLisandro Dalcin jac->maxrowsum = .9; 76316d9e3a6SLisandro Dalcin jac->coarsentype = 6; 76416d9e3a6SLisandro Dalcin jac->measuretype = 0; 7650f1074feSSatish Balay jac->gridsweeps[0] = jac->gridsweeps[1] = jac->gridsweeps[2] = 1; 7668f87f92bSBarry Smith jac->relaxtype[0] = jac->relaxtype[1] = 6; /* Defaults to SYMMETRIC since in PETSc we are using a a PC - most likely with CG */ 7670f1074feSSatish Balay jac->relaxtype[2] = 9; /*G.E. */ 76816d9e3a6SLisandro Dalcin jac->relaxweight = 1.0; 76916d9e3a6SLisandro Dalcin jac->outerrelaxweight = 1.0; 77016d9e3a6SLisandro Dalcin jac->relaxorder = 1; 7710f1074feSSatish Balay jac->interptype = 0; 7720f1074feSSatish Balay jac->agg_nl = 0; 7730f1074feSSatish Balay jac->pmax = 0; 7740f1074feSSatish Balay jac->truncfactor = 0.0; 7750f1074feSSatish Balay jac->agg_num_paths = 1; 7768f87f92bSBarry Smith 7778f87f92bSBarry Smith jac->nodal_coarsen = 0; 7788f87f92bSBarry Smith jac->nodal_relax = PETSC_FALSE; 7798f87f92bSBarry Smith jac->nodal_relax_levels = 1; 780fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetCycleType,(jac->hsolver,jac->cycletype)); 781fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetMaxLevels,(jac->hsolver,jac->maxlevels)); 782fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetMaxIter,(jac->hsolver,jac->maxiter)); 783fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetTol,(jac->hsolver,jac->tol)); 784fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetTruncFactor,(jac->hsolver,jac->truncfactor)); 785fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetStrongThreshold,(jac->hsolver,jac->strongthreshold)); 786fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetMaxRowSum,(jac->hsolver,jac->maxrowsum)); 787fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetCoarsenType,(jac->hsolver,jac->coarsentype)); 788fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetMeasureType,(jac->hsolver,jac->measuretype)); 789fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxOrder,(jac->hsolver, jac->relaxorder)); 790fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetInterpType,(jac->hsolver,jac->interptype)); 791fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetAggNumLevels,(jac->hsolver,jac->agg_nl)); 792fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetPMaxElmts,(jac->hsolver,jac->pmax)); 793fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetNumPaths,(jac->hsolver,jac->agg_num_paths)); 794fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetRelaxType,(jac->hsolver, jac->relaxtype[0])); /*defaults coarse to 9*/ 795fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_BoomerAMGSetNumSweeps,(jac->hsolver, jac->gridsweeps[0])); /*defaults coarse to 1 */ 79616d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 79716d9e3a6SLisandro Dalcin } 798503cfb0cSBarry Smith ierr = PetscFree(jac->hypre_type);CHKERRQ(ierr); 7992fa5cd67SKarl Rupp 8000298fd71SBarry Smith jac->hypre_type = NULL; 801*22b6d1caSBarry Smith SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown HYPRE preconditioner %s; Choices are pilut, parasails, boomeramg",name); 80216d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 80316d9e3a6SLisandro Dalcin } 80416d9e3a6SLisandro Dalcin 80516d9e3a6SLisandro Dalcin /* 80616d9e3a6SLisandro Dalcin It only gets here if the HYPRE type has not been set before the call to 80716d9e3a6SLisandro Dalcin ...SetFromOptions() which actually is most of the time 80816d9e3a6SLisandro Dalcin */ 80916d9e3a6SLisandro Dalcin #undef __FUNCT__ 81016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCSetFromOptions_HYPRE" 811e55864a3SBarry Smith static PetscErrorCode PCSetFromOptions_HYPRE(PetscOptionsObjectType *PetscOptionsObject,PC pc) 81216d9e3a6SLisandro Dalcin { 81316d9e3a6SLisandro Dalcin PetscErrorCode ierr; 8144ddd07fcSJed Brown PetscInt indx; 815*22b6d1caSBarry Smith const char *type[] = {"pilut","parasails","boomeramg"}; 816ace3abfcSBarry Smith PetscBool flg; 81716d9e3a6SLisandro Dalcin 81816d9e3a6SLisandro Dalcin PetscFunctionBegin; 819e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HYPRE preconditioner options");CHKERRQ(ierr); 820*22b6d1caSBarry Smith ierr = PetscOptionsEList("-pc_hypre_type","HYPRE preconditioner type","PCHYPRESetType",type,3,"boomeramg",&indx,&flg);CHKERRQ(ierr); 82116d9e3a6SLisandro Dalcin if (flg) { 82216d9e3a6SLisandro Dalcin ierr = PCHYPRESetType_HYPRE(pc,type[indx]);CHKERRQ(ierr); 82302a17cd4SBarry Smith } else { 82402a17cd4SBarry Smith ierr = PCHYPRESetType_HYPRE(pc,"boomeramg");CHKERRQ(ierr); 82516d9e3a6SLisandro Dalcin } 82616d9e3a6SLisandro Dalcin if (pc->ops->setfromoptions) { 8273931853cSBarry Smith ierr = pc->ops->setfromoptions(PetscOptionsObject,pc);CHKERRQ(ierr); 82816d9e3a6SLisandro Dalcin } 82916d9e3a6SLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 83016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 83116d9e3a6SLisandro Dalcin } 83216d9e3a6SLisandro Dalcin 83316d9e3a6SLisandro Dalcin #undef __FUNCT__ 83416d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPRESetType" 83516d9e3a6SLisandro Dalcin /*@C 83616d9e3a6SLisandro Dalcin PCHYPRESetType - Sets which hypre preconditioner you wish to use 83716d9e3a6SLisandro Dalcin 83816d9e3a6SLisandro Dalcin Input Parameters: 83916d9e3a6SLisandro Dalcin + pc - the preconditioner context 840*22b6d1caSBarry Smith - name - either pilut, parasails, boomeramg 84116d9e3a6SLisandro Dalcin 84216d9e3a6SLisandro Dalcin Options Database Keys: 843*22b6d1caSBarry Smith -pc_hypre_type - One of pilut, parasails, boomeramg 84416d9e3a6SLisandro Dalcin 84516d9e3a6SLisandro Dalcin Level: intermediate 84616d9e3a6SLisandro Dalcin 84716d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 84816d9e3a6SLisandro Dalcin PCHYPRE 84916d9e3a6SLisandro Dalcin 85016d9e3a6SLisandro Dalcin @*/ 8517087cfbeSBarry Smith PetscErrorCode PCHYPRESetType(PC pc,const char name[]) 85216d9e3a6SLisandro Dalcin { 8534ac538c5SBarry Smith PetscErrorCode ierr; 85416d9e3a6SLisandro Dalcin 85516d9e3a6SLisandro Dalcin PetscFunctionBegin; 8560700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 85716d9e3a6SLisandro Dalcin PetscValidCharPointer(name,2); 8584ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCHYPRESetType_C",(PC,const char[]),(pc,name));CHKERRQ(ierr); 85916d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 86016d9e3a6SLisandro Dalcin } 86116d9e3a6SLisandro Dalcin 86216d9e3a6SLisandro Dalcin #undef __FUNCT__ 86316d9e3a6SLisandro Dalcin #define __FUNCT__ "PCHYPREGetType" 86416d9e3a6SLisandro Dalcin /*@C 86516d9e3a6SLisandro Dalcin PCHYPREGetType - Gets which hypre preconditioner you are using 86616d9e3a6SLisandro Dalcin 86716d9e3a6SLisandro Dalcin Input Parameter: 86816d9e3a6SLisandro Dalcin . pc - the preconditioner context 86916d9e3a6SLisandro Dalcin 87016d9e3a6SLisandro Dalcin Output Parameter: 871*22b6d1caSBarry Smith . name - either pilut, parasails, boomeramg 87216d9e3a6SLisandro Dalcin 87316d9e3a6SLisandro Dalcin Level: intermediate 87416d9e3a6SLisandro Dalcin 87516d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCHYPRESetType(), PCType (for list of available types), PC, 87616d9e3a6SLisandro Dalcin PCHYPRE 87716d9e3a6SLisandro Dalcin 87816d9e3a6SLisandro Dalcin @*/ 8797087cfbeSBarry Smith PetscErrorCode PCHYPREGetType(PC pc,const char *name[]) 88016d9e3a6SLisandro Dalcin { 8814ac538c5SBarry Smith PetscErrorCode ierr; 88216d9e3a6SLisandro Dalcin 88316d9e3a6SLisandro Dalcin PetscFunctionBegin; 8840700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 88516d9e3a6SLisandro Dalcin PetscValidPointer(name,2); 8864ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCHYPREGetType_C",(PC,const char*[]),(pc,name));CHKERRQ(ierr); 88716d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 88816d9e3a6SLisandro Dalcin } 88916d9e3a6SLisandro Dalcin 89016d9e3a6SLisandro Dalcin /*MC 89116d9e3a6SLisandro Dalcin PCHYPRE - Allows you to use the matrix element based preconditioners in the LLNL package hypre 89216d9e3a6SLisandro Dalcin 89316d9e3a6SLisandro Dalcin Options Database Keys: 894*22b6d1caSBarry Smith + -pc_hypre_type - One of pilut, parasails, boomeramg 89516d9e3a6SLisandro Dalcin - Too many others to list, run with -pc_type hypre -pc_hypre_type XXX -help to see options for the XXX 89616d9e3a6SLisandro Dalcin preconditioner 89716d9e3a6SLisandro Dalcin 89816d9e3a6SLisandro Dalcin Level: intermediate 89916d9e3a6SLisandro Dalcin 90016d9e3a6SLisandro Dalcin Notes: Apart from pc_hypre_type (for which there is PCHYPRESetType()), 90116d9e3a6SLisandro Dalcin the many hypre options can ONLY be set via the options database (e.g. the command line 90216d9e3a6SLisandro Dalcin or with PetscOptionsSetValue(), there are no functions to set them) 90316d9e3a6SLisandro Dalcin 90416d9e3a6SLisandro Dalcin The options -pc_hypre_boomeramg_max_iter and -pc_hypre_boomeramg_rtol refer to the number of iterations 9050f1074feSSatish Balay (V-cycles) and tolerance that boomeramg does EACH time it is called. So for example, if 9060f1074feSSatish Balay -pc_hypre_boomeramg_max_iter is set to 2 then 2-V-cycles are being used to define the preconditioner 9070f1074feSSatish Balay (-pc_hypre_boomeramg_rtol should be set to 0.0 - the default - to strictly use a fixed number of 9088f87f92bSBarry Smith iterations per hypre call). -ksp_max_it and -ksp_rtol STILL determine the total number of iterations 9090f1074feSSatish Balay and tolerance for the Krylov solver. For example, if -pc_hypre_boomeramg_max_iter is 2 and -ksp_max_it is 10 9100f1074feSSatish Balay then AT MOST twenty V-cycles of boomeramg will be called. 91116d9e3a6SLisandro Dalcin 9120f1074feSSatish Balay Note that the option -pc_hypre_boomeramg_relax_type_all defaults to symmetric relaxation 9130f1074feSSatish Balay (symmetric-SOR/Jacobi), which is required for Krylov solvers like CG that expect symmetry. 9140f1074feSSatish Balay Otherwise, you may want to use -pc_hypre_boomeramg_relax_type_all SOR/Jacobi. 91516d9e3a6SLisandro Dalcin If you wish to use BoomerAMG WITHOUT a Krylov method use -ksp_type richardson NOT -ksp_type preonly 91616d9e3a6SLisandro Dalcin and use -ksp_max_it to control the number of V-cycles. 91716d9e3a6SLisandro Dalcin (see the PETSc FAQ.html at the PETSc website under the Documentation tab). 91816d9e3a6SLisandro Dalcin 91916d9e3a6SLisandro Dalcin 2007-02-03 Using HYPRE-1.11.1b, the routine HYPRE_BoomerAMGSolveT and the option 92016d9e3a6SLisandro Dalcin -pc_hypre_parasails_reuse were failing with SIGSEGV. Dalcin L. 92116d9e3a6SLisandro Dalcin 9229e5bc791SBarry Smith See PCPFMG for access to the hypre Struct PFMG solver 9239e5bc791SBarry Smith 92416d9e3a6SLisandro Dalcin .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 9259e5bc791SBarry Smith PCHYPRESetType(), PCPFMG 92616d9e3a6SLisandro Dalcin 92716d9e3a6SLisandro Dalcin M*/ 92816d9e3a6SLisandro Dalcin 92916d9e3a6SLisandro Dalcin #undef __FUNCT__ 93016d9e3a6SLisandro Dalcin #define __FUNCT__ "PCCreate_HYPRE" 9318cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_HYPRE(PC pc) 93216d9e3a6SLisandro Dalcin { 93316d9e3a6SLisandro Dalcin PC_HYPRE *jac; 93416d9e3a6SLisandro Dalcin PetscErrorCode ierr; 93516d9e3a6SLisandro Dalcin 93616d9e3a6SLisandro Dalcin PetscFunctionBegin; 937b00a9115SJed Brown ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr); 9382fa5cd67SKarl Rupp 93916d9e3a6SLisandro Dalcin pc->data = jac; 94016d9e3a6SLisandro Dalcin pc->ops->destroy = PCDestroy_HYPRE; 94116d9e3a6SLisandro Dalcin pc->ops->setfromoptions = PCSetFromOptions_HYPRE; 94216d9e3a6SLisandro Dalcin pc->ops->setup = PCSetUp_HYPRE; 94316d9e3a6SLisandro Dalcin pc->ops->apply = PCApply_HYPRE; 94416d9e3a6SLisandro Dalcin jac->comm_hypre = MPI_COMM_NULL; 9450298fd71SBarry Smith jac->hypre_type = NULL; 94616d9e3a6SLisandro Dalcin /* duplicate communicator for hypre */ 947ce94432eSBarry Smith ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(jac->comm_hypre));CHKERRQ(ierr); 948bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHYPRESetType_C",PCHYPRESetType_HYPRE);CHKERRQ(ierr); 949bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHYPREGetType_C",PCHYPREGetType_HYPRE);CHKERRQ(ierr); 95016d9e3a6SLisandro Dalcin PetscFunctionReturn(0); 95116d9e3a6SLisandro Dalcin } 952ebc551c0SBarry Smith 953f91d8e95SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------*/ 954f91d8e95SBarry Smith 955b862ddfaSBarry Smith /* this include is needed ONLY to allow access to the private data inside the Mat object specific to hypre */ 956b45d2f2cSJed Brown #include <petsc-private/matimpl.h> 957ebc551c0SBarry Smith 958ebc551c0SBarry Smith typedef struct { 95968326731SBarry Smith MPI_Comm hcomm; /* does not share comm with HYPRE_StructMatrix because need to create solver before getting matrix */ 960f91d8e95SBarry Smith HYPRE_StructSolver hsolver; 9619e5bc791SBarry Smith 9629e5bc791SBarry Smith /* keep copy of PFMG options used so may view them */ 9634ddd07fcSJed Brown PetscInt its; 9649e5bc791SBarry Smith double tol; 9654ddd07fcSJed Brown PetscInt relax_type; 9664ddd07fcSJed Brown PetscInt rap_type; 9674ddd07fcSJed Brown PetscInt num_pre_relax,num_post_relax; 9684ddd07fcSJed Brown PetscInt max_levels; 969ebc551c0SBarry Smith } PC_PFMG; 970ebc551c0SBarry Smith 971ebc551c0SBarry Smith #undef __FUNCT__ 972ebc551c0SBarry Smith #define __FUNCT__ "PCDestroy_PFMG" 973ebc551c0SBarry Smith PetscErrorCode PCDestroy_PFMG(PC pc) 974ebc551c0SBarry Smith { 975ebc551c0SBarry Smith PetscErrorCode ierr; 976f91d8e95SBarry Smith PC_PFMG *ex = (PC_PFMG*) pc->data; 977ebc551c0SBarry Smith 978ebc551c0SBarry Smith PetscFunctionBegin; 9792fa5cd67SKarl Rupp if (ex->hsolver) PetscStackCallStandard(HYPRE_StructPFMGDestroy,(ex->hsolver)); 980f91d8e95SBarry Smith ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr); 981c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 982ebc551c0SBarry Smith PetscFunctionReturn(0); 983ebc551c0SBarry Smith } 984ebc551c0SBarry Smith 9859e5bc791SBarry Smith static const char *PFMGRelaxType[] = {"Jacobi","Weighted-Jacobi","symmetric-Red/Black-Gauss-Seidel","Red/Black-Gauss-Seidel"}; 9869e5bc791SBarry Smith static const char *PFMGRAPType[] = {"Galerkin","non-Galerkin"}; 9879e5bc791SBarry Smith 988ebc551c0SBarry Smith #undef __FUNCT__ 989ebc551c0SBarry Smith #define __FUNCT__ "PCView_PFMG" 990ebc551c0SBarry Smith PetscErrorCode PCView_PFMG(PC pc,PetscViewer viewer) 991ebc551c0SBarry Smith { 992ebc551c0SBarry Smith PetscErrorCode ierr; 993ace3abfcSBarry Smith PetscBool iascii; 994f91d8e95SBarry Smith PC_PFMG *ex = (PC_PFMG*) pc->data; 995ebc551c0SBarry Smith 996ebc551c0SBarry Smith PetscFunctionBegin; 997251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 9989e5bc791SBarry Smith if (iascii) { 9999e5bc791SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE PFMG preconditioning\n");CHKERRQ(ierr); 10009e5bc791SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE PFMG: max iterations %d\n",ex->its);CHKERRQ(ierr); 10019e5bc791SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE PFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr); 10029e5bc791SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE PFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr); 10039e5bc791SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE PFMG: RAP type %s\n",PFMGRAPType[ex->rap_type]);CHKERRQ(ierr); 10049e5bc791SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," HYPRE PFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr); 10053b46a515SGlenn Hammond ierr = PetscViewerASCIIPrintf(viewer," HYPRE PFMG: max levels %d\n",ex->max_levels);CHKERRQ(ierr); 10069e5bc791SBarry Smith } 1007ebc551c0SBarry Smith PetscFunctionReturn(0); 1008ebc551c0SBarry Smith } 1009ebc551c0SBarry Smith 10109e5bc791SBarry Smith 1011ebc551c0SBarry Smith #undef __FUNCT__ 1012ebc551c0SBarry Smith #define __FUNCT__ "PCSetFromOptions_PFMG" 1013e55864a3SBarry Smith PetscErrorCode PCSetFromOptions_PFMG(PetscOptionsObjectType *PetscOptionsObject,PC pc) 1014ebc551c0SBarry Smith { 1015ebc551c0SBarry Smith PetscErrorCode ierr; 1016f91d8e95SBarry Smith PC_PFMG *ex = (PC_PFMG*) pc->data; 1017ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 1018ebc551c0SBarry Smith 1019ebc551c0SBarry Smith PetscFunctionBegin; 1020e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"PFMG options");CHKERRQ(ierr); 10210298fd71SBarry Smith ierr = PetscOptionsBool("-pc_pfmg_print_statistics","Print statistics","HYPRE_StructPFMGSetPrintLevel",flg,&flg,NULL);CHKERRQ(ierr); 102268326731SBarry Smith if (flg) { 1023a0324ebeSBarry Smith int level=3; 1024fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSetPrintLevel,(ex->hsolver,level)); 102568326731SBarry Smith } 10260298fd71SBarry Smith ierr = PetscOptionsInt("-pc_pfmg_its","Number of iterations of PFMG to use as preconditioner","HYPRE_StructPFMGSetMaxIter",ex->its,&ex->its,NULL);CHKERRQ(ierr); 1027fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(ex->hsolver,ex->its)); 10280298fd71SBarry Smith ierr = PetscOptionsInt("-pc_pfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_StructPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,NULL);CHKERRQ(ierr); 1029fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSetNumPreRelax,(ex->hsolver,ex->num_pre_relax)); 10300298fd71SBarry Smith ierr = PetscOptionsInt("-pc_pfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_StructPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,NULL);CHKERRQ(ierr); 1031fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSetNumPostRelax,(ex->hsolver,ex->num_post_relax)); 10329e5bc791SBarry Smith 10330298fd71SBarry Smith ierr = PetscOptionsInt("-pc_pfmg_max_levels","Max Levels for MG hierarchy","HYPRE_StructPFMGSetMaxLevels",ex->max_levels,&ex->max_levels,NULL);CHKERRQ(ierr); 1034fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSetMaxLevels,(ex->hsolver,ex->max_levels)); 10353b46a515SGlenn Hammond 10360298fd71SBarry Smith ierr = PetscOptionsReal("-pc_pfmg_tol","Tolerance of PFMG","HYPRE_StructPFMGSetTol",ex->tol,&ex->tol,NULL);CHKERRQ(ierr); 1037fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSetTol,(ex->hsolver,ex->tol)); 10380298fd71SBarry Smith ierr = PetscOptionsEList("-pc_pfmg_relax_type","Relax type for the up and down cycles","HYPRE_StructPFMGSetRelaxType",PFMGRelaxType,ALEN(PFMGRelaxType),PFMGRelaxType[ex->relax_type],&ex->relax_type,NULL);CHKERRQ(ierr); 1039fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSetRelaxType,(ex->hsolver, ex->relax_type)); 10400298fd71SBarry Smith ierr = PetscOptionsEList("-pc_pfmg_rap_type","RAP type","HYPRE_StructPFMGSetRAPType",PFMGRAPType,ALEN(PFMGRAPType),PFMGRAPType[ex->rap_type],&ex->rap_type,NULL);CHKERRQ(ierr); 1041fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSetRAPType,(ex->hsolver, ex->rap_type)); 1042ebc551c0SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 1043ebc551c0SBarry Smith PetscFunctionReturn(0); 1044ebc551c0SBarry Smith } 1045ebc551c0SBarry Smith 1046f91d8e95SBarry Smith #undef __FUNCT__ 1047f91d8e95SBarry Smith #define __FUNCT__ "PCApply_PFMG" 1048f91d8e95SBarry Smith PetscErrorCode PCApply_PFMG(PC pc,Vec x,Vec y) 1049f91d8e95SBarry Smith { 1050f91d8e95SBarry Smith PetscErrorCode ierr; 1051f91d8e95SBarry Smith PC_PFMG *ex = (PC_PFMG*) pc->data; 1052f91d8e95SBarry Smith PetscScalar *xx,*yy; 10534ddd07fcSJed Brown PetscInt ilower[3],iupper[3]; 105468326731SBarry Smith Mat_HYPREStruct *mx = (Mat_HYPREStruct*)(pc->pmat->data); 1055f91d8e95SBarry Smith 1056f91d8e95SBarry Smith PetscFunctionBegin; 1057dff31646SBarry Smith ierr = PetscCitationsRegister(hypreCitation,&cite);CHKERRQ(ierr); 1058aa219208SBarry Smith ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr); 1059f91d8e95SBarry Smith iupper[0] += ilower[0] - 1; 1060f91d8e95SBarry Smith iupper[1] += ilower[1] - 1; 1061f91d8e95SBarry Smith iupper[2] += ilower[2] - 1; 1062f91d8e95SBarry Smith 1063f91d8e95SBarry Smith /* copy x values over to hypre */ 1064fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructVectorSetConstantValues,(mx->hb,0.0)); 1065f91d8e95SBarry Smith ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 10668b1f7689SBarry Smith PetscStackCallStandard(HYPRE_StructVectorSetBoxValues,(mx->hb,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,xx)); 1067f91d8e95SBarry Smith ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 1068fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructVectorAssemble,(mx->hb)); 1069fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSolve,(ex->hsolver,mx->hmat,mx->hb,mx->hx)); 1070f91d8e95SBarry Smith 1071f91d8e95SBarry Smith /* copy solution values back to PETSc */ 1072f91d8e95SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 10738b1f7689SBarry Smith PetscStackCallStandard(HYPRE_StructVectorGetBoxValues,(mx->hx,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,yy)); 1074f91d8e95SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 1075f91d8e95SBarry Smith PetscFunctionReturn(0); 1076f91d8e95SBarry Smith } 1077f91d8e95SBarry Smith 10789e5bc791SBarry Smith #undef __FUNCT__ 10799e5bc791SBarry Smith #define __FUNCT__ "PCApplyRichardson_PFMG" 1080ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_PFMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason) 10819e5bc791SBarry Smith { 10829e5bc791SBarry Smith PC_PFMG *jac = (PC_PFMG*)pc->data; 10839e5bc791SBarry Smith PetscErrorCode ierr; 10844ddd07fcSJed Brown PetscInt oits; 10859e5bc791SBarry Smith 10869e5bc791SBarry Smith PetscFunctionBegin; 1087dff31646SBarry Smith ierr = PetscCitationsRegister(hypreCitation,&cite);CHKERRQ(ierr); 1088fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(jac->hsolver,its*jac->its)); 1089fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSetTol,(jac->hsolver,rtol)); 10909e5bc791SBarry Smith 10919e5bc791SBarry Smith ierr = PCApply_PFMG(pc,b,y);CHKERRQ(ierr); 10928b1f7689SBarry Smith PetscStackCallStandard(HYPRE_StructPFMGGetNumIterations,(jac->hsolver,(HYPRE_Int *)&oits)); 10939e5bc791SBarry Smith *outits = oits; 10949e5bc791SBarry Smith if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS; 10959e5bc791SBarry Smith else *reason = PCRICHARDSON_CONVERGED_RTOL; 1096fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSetTol,(jac->hsolver,jac->tol)); 1097fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSetMaxIter,(jac->hsolver,jac->its)); 10989e5bc791SBarry Smith PetscFunctionReturn(0); 10999e5bc791SBarry Smith } 11009e5bc791SBarry Smith 11019e5bc791SBarry Smith 11023a32d3dbSGlenn Hammond #undef __FUNCT__ 11033a32d3dbSGlenn Hammond #define __FUNCT__ "PCSetUp_PFMG" 11043a32d3dbSGlenn Hammond PetscErrorCode PCSetUp_PFMG(PC pc) 11053a32d3dbSGlenn Hammond { 11063a32d3dbSGlenn Hammond PetscErrorCode ierr; 11073a32d3dbSGlenn Hammond PC_PFMG *ex = (PC_PFMG*) pc->data; 11083a32d3dbSGlenn Hammond Mat_HYPREStruct *mx = (Mat_HYPREStruct*)(pc->pmat->data); 1109ace3abfcSBarry Smith PetscBool flg; 11103a32d3dbSGlenn Hammond 11113a32d3dbSGlenn Hammond PetscFunctionBegin; 1112251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESTRUCT,&flg);CHKERRQ(ierr); 1113ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESTRUCT with this preconditioner"); 11143a32d3dbSGlenn Hammond 11153a32d3dbSGlenn Hammond /* create the hypre solver object and set its information */ 11162fa5cd67SKarl Rupp if (ex->hsolver) PetscStackCallStandard(HYPRE_StructPFMGDestroy,(ex->hsolver)); 1117fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver)); 1118fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSetup,(ex->hsolver,mx->hmat,mx->hb,mx->hx)); 1119fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGSetZeroGuess,(ex->hsolver)); 11203a32d3dbSGlenn Hammond PetscFunctionReturn(0); 11213a32d3dbSGlenn Hammond } 11223a32d3dbSGlenn Hammond 1123ebc551c0SBarry Smith 1124ebc551c0SBarry Smith /*MC 1125ebc551c0SBarry Smith PCPFMG - the hypre PFMG multigrid solver 1126ebc551c0SBarry Smith 1127ebc551c0SBarry Smith Level: advanced 1128ebc551c0SBarry Smith 11299e5bc791SBarry Smith Options Database: 11309e5bc791SBarry Smith + -pc_pfmg_its <its> number of iterations of PFMG to use as preconditioner 11319e5bc791SBarry Smith . -pc_pfmg_num_pre_relax <steps> number of smoothing steps before coarse grid 11329e5bc791SBarry Smith . -pc_pfmg_num_post_relax <steps> number of smoothing steps after coarse grid 11339e5bc791SBarry Smith . -pc_pfmg_tol <tol> tolerance of PFMG 11349e5bc791SBarry Smith . -pc_pfmg_relax_type -relaxation type for the up and down cycles, one of Jacobi,Weighted-Jacobi,symmetric-Red/Black-Gauss-Seidel,Red/Black-Gauss-Seidel 11359e5bc791SBarry Smith - -pc_pfmg_rap_type - type of coarse matrix generation, one of Galerkin,non-Galerkin 1136f91d8e95SBarry Smith 11379e5bc791SBarry Smith Notes: This is for CELL-centered descretizations 11389e5bc791SBarry Smith 11398e395302SJed Brown This must be used with the MATHYPRESTRUCT matrix type. 1140aa219208SBarry Smith This is less general than in hypre, it supports only one block per process defined by a PETSc DMDA. 11419e5bc791SBarry Smith 11429e5bc791SBarry Smith .seealso: PCMG, MATHYPRESTRUCT 1143ebc551c0SBarry Smith M*/ 1144ebc551c0SBarry Smith 1145ebc551c0SBarry Smith #undef __FUNCT__ 1146ebc551c0SBarry Smith #define __FUNCT__ "PCCreate_PFMG" 11478cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_PFMG(PC pc) 1148ebc551c0SBarry Smith { 1149ebc551c0SBarry Smith PetscErrorCode ierr; 1150ebc551c0SBarry Smith PC_PFMG *ex; 1151ebc551c0SBarry Smith 1152ebc551c0SBarry Smith PetscFunctionBegin; 1153b00a9115SJed Brown ierr = PetscNew(&ex);CHKERRQ(ierr); \ 115468326731SBarry Smith pc->data = ex; 1155ebc551c0SBarry Smith 11569e5bc791SBarry Smith ex->its = 1; 11579e5bc791SBarry Smith ex->tol = 1.e-8; 11589e5bc791SBarry Smith ex->relax_type = 1; 11599e5bc791SBarry Smith ex->rap_type = 0; 11609e5bc791SBarry Smith ex->num_pre_relax = 1; 11619e5bc791SBarry Smith ex->num_post_relax = 1; 11623b46a515SGlenn Hammond ex->max_levels = 0; 11639e5bc791SBarry Smith 1164ebc551c0SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_PFMG; 1165ebc551c0SBarry Smith pc->ops->view = PCView_PFMG; 1166ebc551c0SBarry Smith pc->ops->destroy = PCDestroy_PFMG; 1167f91d8e95SBarry Smith pc->ops->apply = PCApply_PFMG; 11689e5bc791SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_PFMG; 116968326731SBarry Smith pc->ops->setup = PCSetUp_PFMG; 11702fa5cd67SKarl Rupp 1171ce94432eSBarry Smith ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ex->hcomm));CHKERRQ(ierr); 1172fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_StructPFMGCreate,(ex->hcomm,&ex->hsolver)); 1173ebc551c0SBarry Smith PetscFunctionReturn(0); 1174ebc551c0SBarry Smith } 1175d851a50bSGlenn Hammond 1176325fc9f4SBarry Smith /* ---------------------------------------------------------------------------------------------------------------------------------------------------*/ 1177325fc9f4SBarry Smith 1178d851a50bSGlenn Hammond /* we know we are working with a HYPRE_SStructMatrix */ 1179d851a50bSGlenn Hammond typedef struct { 1180d851a50bSGlenn Hammond MPI_Comm hcomm; /* does not share comm with HYPRE_SStructMatrix because need to create solver before getting matrix */ 1181d851a50bSGlenn Hammond HYPRE_SStructSolver ss_solver; 1182d851a50bSGlenn Hammond 1183d851a50bSGlenn Hammond /* keep copy of SYSPFMG options used so may view them */ 11844ddd07fcSJed Brown PetscInt its; 1185d851a50bSGlenn Hammond double tol; 11864ddd07fcSJed Brown PetscInt relax_type; 11874ddd07fcSJed Brown PetscInt num_pre_relax,num_post_relax; 1188d851a50bSGlenn Hammond } PC_SysPFMG; 1189d851a50bSGlenn Hammond 1190d851a50bSGlenn Hammond #undef __FUNCT__ 1191d851a50bSGlenn Hammond #define __FUNCT__ "PCDestroy_SysPFMG" 1192d851a50bSGlenn Hammond PetscErrorCode PCDestroy_SysPFMG(PC pc) 1193d851a50bSGlenn Hammond { 1194d851a50bSGlenn Hammond PetscErrorCode ierr; 1195d851a50bSGlenn Hammond PC_SysPFMG *ex = (PC_SysPFMG*) pc->data; 1196d851a50bSGlenn Hammond 1197d851a50bSGlenn Hammond PetscFunctionBegin; 11982fa5cd67SKarl Rupp if (ex->ss_solver) PetscStackCallStandard(HYPRE_SStructSysPFMGDestroy,(ex->ss_solver)); 1199d851a50bSGlenn Hammond ierr = MPI_Comm_free(&ex->hcomm);CHKERRQ(ierr); 1200c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1201d851a50bSGlenn Hammond PetscFunctionReturn(0); 1202d851a50bSGlenn Hammond } 1203d851a50bSGlenn Hammond 1204d851a50bSGlenn Hammond static const char *SysPFMGRelaxType[] = {"Weighted-Jacobi","Red/Black-Gauss-Seidel"}; 1205d851a50bSGlenn Hammond 1206d851a50bSGlenn Hammond #undef __FUNCT__ 1207d851a50bSGlenn Hammond #define __FUNCT__ "PCView_SysPFMG" 1208d851a50bSGlenn Hammond PetscErrorCode PCView_SysPFMG(PC pc,PetscViewer viewer) 1209d851a50bSGlenn Hammond { 1210d851a50bSGlenn Hammond PetscErrorCode ierr; 1211ace3abfcSBarry Smith PetscBool iascii; 1212d851a50bSGlenn Hammond PC_SysPFMG *ex = (PC_SysPFMG*) pc->data; 1213d851a50bSGlenn Hammond 1214d851a50bSGlenn Hammond PetscFunctionBegin; 1215251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1216d851a50bSGlenn Hammond if (iascii) { 1217d851a50bSGlenn Hammond ierr = PetscViewerASCIIPrintf(viewer," HYPRE SysPFMG preconditioning\n");CHKERRQ(ierr); 1218d851a50bSGlenn Hammond ierr = PetscViewerASCIIPrintf(viewer," HYPRE SysPFMG: max iterations %d\n",ex->its);CHKERRQ(ierr); 1219d851a50bSGlenn Hammond ierr = PetscViewerASCIIPrintf(viewer," HYPRE SysPFMG: tolerance %g\n",ex->tol);CHKERRQ(ierr); 1220d851a50bSGlenn Hammond ierr = PetscViewerASCIIPrintf(viewer," HYPRE SysPFMG: relax type %s\n",PFMGRelaxType[ex->relax_type]);CHKERRQ(ierr); 1221d851a50bSGlenn Hammond ierr = PetscViewerASCIIPrintf(viewer," HYPRE SysPFMG: number pre-relax %d post-relax %d\n",ex->num_pre_relax,ex->num_post_relax);CHKERRQ(ierr); 1222d851a50bSGlenn Hammond } 1223d851a50bSGlenn Hammond PetscFunctionReturn(0); 1224d851a50bSGlenn Hammond } 1225d851a50bSGlenn Hammond 1226d851a50bSGlenn Hammond 1227d851a50bSGlenn Hammond #undef __FUNCT__ 1228d851a50bSGlenn Hammond #define __FUNCT__ "PCSetFromOptions_SysPFMG" 1229e55864a3SBarry Smith PetscErrorCode PCSetFromOptions_SysPFMG(PetscOptionsObjectType *PetscOptionsObject,PC pc) 1230d851a50bSGlenn Hammond { 1231d851a50bSGlenn Hammond PetscErrorCode ierr; 1232d851a50bSGlenn Hammond PC_SysPFMG *ex = (PC_SysPFMG*) pc->data; 1233ace3abfcSBarry Smith PetscBool flg = PETSC_FALSE; 1234d851a50bSGlenn Hammond 1235d851a50bSGlenn Hammond PetscFunctionBegin; 1236e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SysPFMG options");CHKERRQ(ierr); 12370298fd71SBarry Smith ierr = PetscOptionsBool("-pc_syspfmg_print_statistics","Print statistics","HYPRE_SStructSysPFMGSetPrintLevel",flg,&flg,NULL);CHKERRQ(ierr); 1238d851a50bSGlenn Hammond if (flg) { 1239d851a50bSGlenn Hammond int level=3; 1240fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGSetPrintLevel,(ex->ss_solver,level)); 1241d851a50bSGlenn Hammond } 12420298fd71SBarry Smith ierr = PetscOptionsInt("-pc_syspfmg_its","Number of iterations of SysPFMG to use as preconditioner","HYPRE_SStructSysPFMGSetMaxIter",ex->its,&ex->its,NULL);CHKERRQ(ierr); 1243fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(ex->ss_solver,ex->its)); 12440298fd71SBarry Smith ierr = PetscOptionsInt("-pc_syspfmg_num_pre_relax","Number of smoothing steps before coarse grid","HYPRE_SStructSysPFMGSetNumPreRelax",ex->num_pre_relax,&ex->num_pre_relax,NULL);CHKERRQ(ierr); 1245fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGSetNumPreRelax,(ex->ss_solver,ex->num_pre_relax)); 12460298fd71SBarry Smith ierr = PetscOptionsInt("-pc_syspfmg_num_post_relax","Number of smoothing steps after coarse grid","HYPRE_SStructSysPFMGSetNumPostRelax",ex->num_post_relax,&ex->num_post_relax,NULL);CHKERRQ(ierr); 1247fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGSetNumPostRelax,(ex->ss_solver,ex->num_post_relax)); 1248d851a50bSGlenn Hammond 12490298fd71SBarry Smith ierr = PetscOptionsReal("-pc_syspfmg_tol","Tolerance of SysPFMG","HYPRE_SStructSysPFMGSetTol",ex->tol,&ex->tol,NULL);CHKERRQ(ierr); 1250fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(ex->ss_solver,ex->tol)); 12510298fd71SBarry Smith ierr = PetscOptionsEList("-pc_syspfmg_relax_type","Relax type for the up and down cycles","HYPRE_SStructSysPFMGSetRelaxType",SysPFMGRelaxType,4,SysPFMGRelaxType[ex->relax_type],&ex->relax_type,NULL);CHKERRQ(ierr); 1252fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGSetRelaxType,(ex->ss_solver, ex->relax_type)); 1253d851a50bSGlenn Hammond ierr = PetscOptionsTail();CHKERRQ(ierr); 1254d851a50bSGlenn Hammond PetscFunctionReturn(0); 1255d851a50bSGlenn Hammond } 1256d851a50bSGlenn Hammond 1257d851a50bSGlenn Hammond #undef __FUNCT__ 1258d851a50bSGlenn Hammond #define __FUNCT__ "PCApply_SysPFMG" 1259d851a50bSGlenn Hammond PetscErrorCode PCApply_SysPFMG(PC pc,Vec x,Vec y) 1260d851a50bSGlenn Hammond { 1261d851a50bSGlenn Hammond PetscErrorCode ierr; 1262d851a50bSGlenn Hammond PC_SysPFMG *ex = (PC_SysPFMG*) pc->data; 1263d851a50bSGlenn Hammond PetscScalar *xx,*yy; 12644ddd07fcSJed Brown PetscInt ilower[3],iupper[3]; 1265d851a50bSGlenn Hammond Mat_HYPRESStruct *mx = (Mat_HYPRESStruct*)(pc->pmat->data); 12664ddd07fcSJed Brown PetscInt ordering= mx->dofs_order; 12674ddd07fcSJed Brown PetscInt nvars = mx->nvars; 12684ddd07fcSJed Brown PetscInt part = 0; 12694ddd07fcSJed Brown PetscInt size; 12704ddd07fcSJed Brown PetscInt i; 1271d851a50bSGlenn Hammond 1272d851a50bSGlenn Hammond PetscFunctionBegin; 1273dff31646SBarry Smith ierr = PetscCitationsRegister(hypreCitation,&cite);CHKERRQ(ierr); 1274aa219208SBarry Smith ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr); 1275d851a50bSGlenn Hammond iupper[0] += ilower[0] - 1; 1276d851a50bSGlenn Hammond iupper[1] += ilower[1] - 1; 1277d851a50bSGlenn Hammond iupper[2] += ilower[2] - 1; 1278d851a50bSGlenn Hammond 1279d851a50bSGlenn Hammond size = 1; 12802fa5cd67SKarl Rupp for (i= 0; i< 3; i++) size *= (iupper[i]-ilower[i]+1); 12812fa5cd67SKarl Rupp 1282d851a50bSGlenn Hammond /* copy x values over to hypre for variable ordering */ 1283d851a50bSGlenn Hammond if (ordering) { 1284fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0)); 1285d851a50bSGlenn Hammond ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 12868b1f7689SBarry Smith for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,xx+(size*i))); 1287d851a50bSGlenn Hammond ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 1288fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b)); 1289fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x)); 1290fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x)); 1291d851a50bSGlenn Hammond 1292d851a50bSGlenn Hammond /* copy solution values back to PETSc */ 1293d851a50bSGlenn Hammond ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 12948b1f7689SBarry Smith for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,yy+(size*i))); 1295d851a50bSGlenn Hammond ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 1296a65764d7SBarry Smith } else { /* nodal ordering must be mapped to variable ordering for sys_pfmg */ 1297d851a50bSGlenn Hammond PetscScalar *z; 12984ddd07fcSJed Brown PetscInt j, k; 1299d851a50bSGlenn Hammond 1300785e854fSJed Brown ierr = PetscMalloc1(nvars*size,&z);CHKERRQ(ierr); 1301fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0)); 1302d851a50bSGlenn Hammond ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 1303d851a50bSGlenn Hammond 1304d851a50bSGlenn Hammond /* transform nodal to hypre's variable ordering for sys_pfmg */ 1305d851a50bSGlenn Hammond for (i= 0; i< size; i++) { 1306d851a50bSGlenn Hammond k= i*nvars; 13072fa5cd67SKarl Rupp for (j= 0; j< nvars; j++) z[j*size+i]= xx[k+j]; 1308d851a50bSGlenn Hammond } 13098b1f7689SBarry Smith for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,z+(size*i))); 1310d851a50bSGlenn Hammond ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 1311fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b)); 1312fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGSolve,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x)); 1313d851a50bSGlenn Hammond 1314d851a50bSGlenn Hammond /* copy solution values back to PETSc */ 1315d851a50bSGlenn Hammond ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 13168b1f7689SBarry Smith for (i= 0; i< nvars; i++) PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,z+(size*i))); 1317d851a50bSGlenn Hammond /* transform hypre's variable ordering for sys_pfmg to nodal ordering */ 1318d851a50bSGlenn Hammond for (i= 0; i< size; i++) { 1319d851a50bSGlenn Hammond k= i*nvars; 13202fa5cd67SKarl Rupp for (j= 0; j< nvars; j++) yy[k+j]= z[j*size+i]; 1321d851a50bSGlenn Hammond } 1322d851a50bSGlenn Hammond ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 1323d851a50bSGlenn Hammond ierr = PetscFree(z);CHKERRQ(ierr); 1324d851a50bSGlenn Hammond } 1325d851a50bSGlenn Hammond PetscFunctionReturn(0); 1326d851a50bSGlenn Hammond } 1327d851a50bSGlenn Hammond 1328d851a50bSGlenn Hammond #undef __FUNCT__ 1329d851a50bSGlenn Hammond #define __FUNCT__ "PCApplyRichardson_SysPFMG" 1330ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_SysPFMG(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason) 1331d851a50bSGlenn Hammond { 1332d851a50bSGlenn Hammond PC_SysPFMG *jac = (PC_SysPFMG*)pc->data; 1333d851a50bSGlenn Hammond PetscErrorCode ierr; 13344ddd07fcSJed Brown PetscInt oits; 1335d851a50bSGlenn Hammond 1336d851a50bSGlenn Hammond PetscFunctionBegin; 1337dff31646SBarry Smith ierr = PetscCitationsRegister(hypreCitation,&cite);CHKERRQ(ierr); 1338fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,its*jac->its)); 1339fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,rtol)); 1340d851a50bSGlenn Hammond ierr = PCApply_SysPFMG(pc,b,y);CHKERRQ(ierr); 13418b1f7689SBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGGetNumIterations,(jac->ss_solver,(HYPRE_Int *)&oits)); 1342d851a50bSGlenn Hammond *outits = oits; 1343d851a50bSGlenn Hammond if (oits == its) *reason = PCRICHARDSON_CONVERGED_ITS; 1344d851a50bSGlenn Hammond else *reason = PCRICHARDSON_CONVERGED_RTOL; 1345fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGSetTol,(jac->ss_solver,jac->tol)); 1346fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGSetMaxIter,(jac->ss_solver,jac->its)); 1347d851a50bSGlenn Hammond PetscFunctionReturn(0); 1348d851a50bSGlenn Hammond } 1349d851a50bSGlenn Hammond 1350d851a50bSGlenn Hammond 1351d851a50bSGlenn Hammond #undef __FUNCT__ 1352d851a50bSGlenn Hammond #define __FUNCT__ "PCSetUp_SysPFMG" 1353d851a50bSGlenn Hammond PetscErrorCode PCSetUp_SysPFMG(PC pc) 1354d851a50bSGlenn Hammond { 1355d851a50bSGlenn Hammond PetscErrorCode ierr; 1356d851a50bSGlenn Hammond PC_SysPFMG *ex = (PC_SysPFMG*) pc->data; 1357d851a50bSGlenn Hammond Mat_HYPRESStruct *mx = (Mat_HYPRESStruct*)(pc->pmat->data); 1358ace3abfcSBarry Smith PetscBool flg; 1359d851a50bSGlenn Hammond 1360d851a50bSGlenn Hammond PetscFunctionBegin; 1361251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATHYPRESSTRUCT,&flg);CHKERRQ(ierr); 1362ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Must use MATHYPRESSTRUCT with this preconditioner"); 1363d851a50bSGlenn Hammond 1364d851a50bSGlenn Hammond /* create the hypre sstruct solver object and set its information */ 13652fa5cd67SKarl Rupp if (ex->ss_solver) PetscStackCallStandard(HYPRE_SStructSysPFMGDestroy,(ex->ss_solver)); 1366fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver)); 1367fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGSetZeroGuess,(ex->ss_solver)); 1368fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGSetup,(ex->ss_solver,mx->ss_mat,mx->ss_b,mx->ss_x)); 1369d851a50bSGlenn Hammond PetscFunctionReturn(0); 1370d851a50bSGlenn Hammond } 1371d851a50bSGlenn Hammond 1372d851a50bSGlenn Hammond 1373d851a50bSGlenn Hammond /*MC 1374d851a50bSGlenn Hammond PCSysPFMG - the hypre SysPFMG multigrid solver 1375d851a50bSGlenn Hammond 1376d851a50bSGlenn Hammond Level: advanced 1377d851a50bSGlenn Hammond 1378d851a50bSGlenn Hammond Options Database: 1379d851a50bSGlenn Hammond + -pc_syspfmg_its <its> number of iterations of SysPFMG to use as preconditioner 1380d851a50bSGlenn Hammond . -pc_syspfmg_num_pre_relax <steps> number of smoothing steps before coarse grid 1381d851a50bSGlenn Hammond . -pc_syspfmg_num_post_relax <steps> number of smoothing steps after coarse grid 1382d851a50bSGlenn Hammond . -pc_syspfmg_tol <tol> tolerance of SysPFMG 1383d851a50bSGlenn Hammond . -pc_syspfmg_relax_type -relaxation type for the up and down cycles, one of Weighted-Jacobi,Red/Black-Gauss-Seidel 1384d851a50bSGlenn Hammond 1385d851a50bSGlenn Hammond Notes: This is for CELL-centered descretizations 1386d851a50bSGlenn Hammond 1387f6680f47SSatish Balay This must be used with the MATHYPRESSTRUCT matrix type. 1388aa219208SBarry Smith This is less general than in hypre, it supports only one part, and one block per process defined by a PETSc DMDA. 1389d851a50bSGlenn Hammond Also, only cell-centered variables. 1390d851a50bSGlenn Hammond 1391d851a50bSGlenn Hammond .seealso: PCMG, MATHYPRESSTRUCT 1392d851a50bSGlenn Hammond M*/ 1393d851a50bSGlenn Hammond 1394d851a50bSGlenn Hammond #undef __FUNCT__ 1395d851a50bSGlenn Hammond #define __FUNCT__ "PCCreate_SysPFMG" 13968cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SysPFMG(PC pc) 1397d851a50bSGlenn Hammond { 1398d851a50bSGlenn Hammond PetscErrorCode ierr; 1399d851a50bSGlenn Hammond PC_SysPFMG *ex; 1400d851a50bSGlenn Hammond 1401d851a50bSGlenn Hammond PetscFunctionBegin; 1402b00a9115SJed Brown ierr = PetscNew(&ex);CHKERRQ(ierr); \ 1403d851a50bSGlenn Hammond pc->data = ex; 1404d851a50bSGlenn Hammond 1405d851a50bSGlenn Hammond ex->its = 1; 1406d851a50bSGlenn Hammond ex->tol = 1.e-8; 1407d851a50bSGlenn Hammond ex->relax_type = 1; 1408d851a50bSGlenn Hammond ex->num_pre_relax = 1; 1409d851a50bSGlenn Hammond ex->num_post_relax = 1; 1410d851a50bSGlenn Hammond 1411d851a50bSGlenn Hammond pc->ops->setfromoptions = PCSetFromOptions_SysPFMG; 1412d851a50bSGlenn Hammond pc->ops->view = PCView_SysPFMG; 1413d851a50bSGlenn Hammond pc->ops->destroy = PCDestroy_SysPFMG; 1414d851a50bSGlenn Hammond pc->ops->apply = PCApply_SysPFMG; 1415d851a50bSGlenn Hammond pc->ops->applyrichardson = PCApplyRichardson_SysPFMG; 1416d851a50bSGlenn Hammond pc->ops->setup = PCSetUp_SysPFMG; 14172fa5cd67SKarl Rupp 1418ce94432eSBarry Smith ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ex->hcomm));CHKERRQ(ierr); 1419fd3f9acdSBarry Smith PetscStackCallStandard(HYPRE_SStructSysPFMGCreate,(ex->hcomm,&ex->ss_solver)); 1420d851a50bSGlenn Hammond PetscFunctionReturn(0); 1421d851a50bSGlenn Hammond } 1422