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